module_PHYSICS_CALLS.F

References to this file elsewhere.
1 !-----------------------------------------------------------------------
2 !
3 !NCEP_MESO:MODEL_LAYER: PHYSICS
4 !
5 !-----------------------------------------------------------------------
6 #include "nmm_loop_basemacros.h"
7 #include "nmm_loop_macros.h"
8 !-----------------------------------------------------------------------
9 !
10       MODULE MODULE_PHYSICS_CALLS
11 !
12 !-----------------------------------------------------------------------
13       USE MODULE_DOMAIN
14       USE MODULE_DM
15       USE MODULE_CONFIGURE
16       USE MODULE_TILES
17       USE MODULE_STATE_DESCRIPTION,ONLY : P_QV,P_QC,P_QR,P_QI,P_QS,P_QG,P_QNI
18       USE MODULE_MODEL_CONSTANTS
19       USE MODULE_RA_GFDLETA,ONLY : CAL_MON_DAY,ZENITH
20       USE MODULE_RADIATION_DRIVER
21       USE MODULE_SF_MYJSFC
22       USE MODULE_SURFACE_DRIVER
23       USE MODULE_PBL_DRIVER
24       USE MODULE_CU_BMJ
25       USE MODULE_CUMULUS_DRIVER
26       USE MODULE_MP_ETANEW
27       USE MODULE_MICROPHYSICS_DRIVER
28       USE MODULE_MICROPHYSICS_ZERO_OUT
29 !-----------------------------------------------------------------------
30 !
31       CONTAINS
32 !
33 !-----------------------------------------------------------------------
34 !***********************************************************************
35       SUBROUTINE RADIATION(NTSD,DT,JULDAY,JULYR,XTIME,JULIAN            &
36      &                    ,IHRST,NPHS,GLAT,GLON                         &
37      &                    ,NRADS,NRADL                                  &
38      &                    ,DETA1,DETA2,AETA1,AETA2,ETA1,ETA2,PDTOP,PT   &
39      &                    ,PD,RES,PINT,T,Q,MOIST,THS,ALBEDO,EPSR        &
40      &                    ,F_ICE,F_RAIN                                 &
41 #ifdef WRF_CHEM
42      &                    ,GD_CLOUD,GD_CLOUD2                           &
43 #endif
44      &                    ,SM,HBM2,LMH,CLDFRA,N_MOIST,RESTRT            &
45      &                    ,RLWTT,RSWTT,RLWIN,RSWIN,RSWINC,RSWOUT        &
46      &                    ,RLWTOA,RSWTOA,CZMEAN                         &
47      &                    ,CFRACL,CFRACM,CFRACH,SIGT4                   &
48      &                    ,ACFRST,NCFRST,ACFRCV,NCFRCV                  &
49      &                    ,CUPPT,VEGFRC,SNOW,HTOP,HBOT                  &
50      &                    ,Z,SICE,NUM_AEROSOLC,NUM_OZMIXM               &
51      &                    ,GRID,CONFIG_FLAGS                            &
52      &                    ,RTHRATEN                                     &
53 #ifdef WRF_CHEM
54      &                    ,PM2_5_DRY, PM2_5_WATER, PM2_5_DRY_EC         &
55      &                    ,TAUAER1, TAUAER2, TAUAER3, TAUAER4           &
56      &                    ,GAER1, GAER2, GAER3, GAER4                   &
57      &                    ,WAER1, WAER2, WAER3, WAER4                   &
58 #endif
59      &                    ,IDS,IDE,JDS,JDE,KDS,KDE                      &
60      &                    ,IMS,IME,JMS,JME,KMS,KME                      &
61      &                    ,ITS,ITE,JTS,JTE,KTS,KTE)
62 !***  NOTE ***
63 ! RLWIN  - downward longwave at the surface (=TOTLWDN, now a local array)
64 ! RSWIN  - downward shortwave at the surface (=TOTSWDN, now a local array)
65 ! RSWINC - CLEAR-SKY downward shortwave at the surface (=TOTSWDNC, new for AQ)
66 !***********************************************************************
67 !$$$  SUBPROGRAM DOCUMENTATION BLOCK
68 !                .      .    .     
69 ! SUBPROGRAM:    RADIATION   RADIATION OUTER DRIVER
70 !   PRGRMMR: BLACK           ORG: W/NP22     DATE: 2002-06-04       
71 !     
72 ! ABSTRACT:
73 !     RADIATION SERVES AS THE INTERFACE BETWEEN THE NCEP NONHYDROSTATIC
74 !     MESOSCALE MODEL AND THE WRF RADIATION DRIVER.
75 !     
76 ! PROGRAM HISTORY LOG:
77 !   02-06-04  BLACK      - ORIGINATOR
78 !   02-09-09  WOLFE      - CONVERTING TO GLOBAL INDEXING
79 !   04-11-18  BLACK      - THREADED
80 !     
81 ! USAGE: CALL RADIATION FROM SOLVE_NMM      
82 !
83 ! ATTRIBUTES:
84 !   LANGUAGE: FORTRAN 90
85 !   MACHINE : IBM 
86 !$$$  
87 !-----------------------------------------------------------------------
88 !
89       IMPLICIT NONE
90 !
91 !-----------------------------------------------------------------------
92 !
93       INTEGER,INTENT(IN) :: IDS,IDE,JDS,JDE,KDS,KDE                     &
94      &                     ,IMS,IME,JMS,JME,KMS,KME                     &
95      &                     ,ITS,ITE,JTS,JTE,KTS,KTE                     &
96      &                     ,IHRST,JULDAY,JULYR                          &
97      &                     ,N_MOIST,NPHS,NRADL,NRADS,NTSD               &
98      &                     ,NUM_AEROSOLC,NUM_OZMIXM
99 !
100       INTEGER,DIMENSION(IMS:IME,JMS:JME),INTENT(INOUT) :: NCFRCV,NCFRST
101 !
102       REAL,INTENT(IN) :: DT,PDTOP,PT,XTIME,JULIAN
103 !
104       INTEGER,DIMENSION(IMS:IME,JMS:JME),INTENT(IN) :: LMH
105 !
106       REAL,DIMENSION(KMS:KME-1),INTENT(IN) :: AETA1,AETA2,DETA1,DETA2
107 !
108       REAL,DIMENSION(KMS:KME),INTENT(IN) :: ETA1,ETA2
109 !
110       REAL,DIMENSION(IMS:IME,JMS:JME),INTENT(IN) :: ALBEDO              &
111      &                                             ,EPSR,GLAT,GLON      &
112      &                                             ,HBM2                &
113      &                                             ,PD,RES,SM           &
114      &                                             ,SNOW,THS,VEGFRC,SICE
115       REAL,DIMENSION(IMS:IME,JMS:JME),INTENT(INOUT) :: CUPPT            
116 
117 !
118       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(IN) :: F_ICE       &
119      &                                                     ,F_RAIN      &
120      &                                                     ,Q,T,Z
121 !
122       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(INOUT) :: RTHRATEN
123 !
124       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME,N_MOIST),INTENT(INOUT) :: MOIST
125 !
126       REAL,DIMENSION(IMS:IME,JMS:JME),INTENT(INOUT) :: ACFRCV,ACFRST    &
127      &                                                ,RLWIN,RLWTOA     &
128      &                                                ,RSWIN,RSWOUT     &
129      &                                                ,HBOT,HTOP        &
130      &                                                ,RSWINC,RSWTOA
131 !
132       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(INOUT) :: PINT     &
133      &                                                        ,RLWTT    &
134      &                                                        ,RSWTT
135 !
136       REAL,DIMENSION(IMS:IME,JMS:JME),INTENT(INOUT) :: CFRACH,CFRACL      &
137      &                                              ,CFRACM,CZMEAN      &
138      &                                              ,SIGT4
139 #ifdef WRF_CHEM
140       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME ),INTENT(IN) ::            &
141      &                              GAER1,GAER2,GAER3,GAER4,            &
142      &                              GD_CLOUD,GD_CLOUD2,                 &
143      &                              PM2_5_DRY,PM2_5_WATER,PM2_5_DRY_EC, &
144      &                              TAUAER1,TAUAER2,TAUAER3,TAUAER4,    &
145      &                              WAER1,WAER2,WAER3,WAER4
146 #endif
147 !
148       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(INOUT) :: CLDFRA
149 !
150       LOGICAL,INTENT(IN) :: RESTRT
151 !
152       TYPE(DOMAIN),TARGET :: GRID
153 !
154       TYPE(GRID_CONFIG_REC_TYPE),INTENT(IN) :: CONFIG_FLAGS
155 !
156 !-----------------------------------------------------------------------
157 !***
158 !***  LOCAL VARIABLES
159 !***
160 !-----------------------------------------------------------------------
161       INTEGER :: I,ICLOUD,IENDX,II,J,JDAY,JMONTH,K,KMNTH,LMHIJ,NRAD
162 !
163       INTEGER,DIMENSION(3) :: IDAT
164       INTEGER,DIMENSION(12) :: MONTH=(/31,28,31,30,31,30,31,31          &
165      &                                ,30,31,30,31/)
166 !
167       REAL :: CAPA,DAYI,DPL,FICE,FRAIN,GMT,HOUR,PDSL,PLYR,PSFC          &
168      &       ,QI,QR,QW,RADT,TIMES,WC,TDUM
169 !
170       REAL,DIMENSION(KMS:KME-1) :: QL,TL
171 !
172       REAL,DIMENSION(IMS:IME,JMS:JME) :: REXNSFC,SWNETDN                &
173      &                                  ,TOT,TSFC,XLAND,XLAT,XLON       &
174      &                                  ,TOTLWDN,TOTSWDN,TOTSWDNC,CZEN  &
175      &                                  ,HBOTR,HTOPR,CUPPTR
176 !
177 !
178       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME) :: DZ,P8W,P_PHY,PI_PHY    &
179      &                                          ,RR,T8W                 &
180      &                                          ,THRATENLW,THRATENSW    &
181      &                                          ,TH_PHY,T_PHY,CLFR
182 !
183 !
184 !***  Different way to include cloud effects in radiation.
185 !
186       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME) :: QC1R,QI1R
187 !
188       LOGICAL :: WARM_RAIN
189 !
190 !-----------------------------------------------------------------------
191 !-----------------------------------------------------------------------
192 !*****
193 !***** NOTE: THIS IS HARDWIRED FOR CALLS TO LONGWAVE AND SHORTWAVE
194 !*****       AT EQUAL INTERVALS
195 !*****
196       NRAD=NRADS
197       RADT=DT*NRADS/60.
198 !-----------------------------------------------------------------------
199 !-----------------------------------------------------------------------
200       CAPA=R_D/CP
201 !-----------------------------------------------------------------------
202 !
203 !$omp parallel do                                                       &
204 !$omp& private(dpl,fice,frain,i,j,k,pdsl,plyr,ql,tl)
205       DO J=MYJS2,MYJE2
206       DO I=MYIS1,MYIE1
207 !
208         PDSL=PD(I,J)*RES(I,J)
209         P8W(I,KTE+1,J)=PT
210         XLAT(I,J)=GLAT(I,J)/DEGRAD
211         XLON(I,J)=GLON(I,J)/DEGRAD
212         XLAND(I,J)=SM(I,J)+1.
213         PSFC=PD(I,J)+PDTOP+PT
214         REXNSFC(I,J)=(PSFC*1.E-5)**CAPA
215         TSFC(I,J)=THS(I,J)*REXNSFC(I,J)
216         T8W(I,1,J)=TSFC(I,J)
217         P8W(I,KTS,J)=ETA1(KTS)*PDTOP+ETA2(KTS)*PDSL+PT
218 !
219 !-----------------------------------------------------------------------
220 !***  FILL THE SINGLE-COLUMN INPUT
221 !-----------------------------------------------------------------------
222 !
223         DO K=KTS,KTE
224           DPL=DETA1(K)*PDTOP+DETA2(K)*PDSL
225           QL(K)=AMAX1(Q(I,K,J),EPSQ)
226           PLYR=AETA1(K)*PDTOP+AETA2(K)*PDSL+PT
227           TL(K)=T(I,K,J)
228 !
229           RR(I,K,J)=PLYR/(R_D*TL(K)*(1.+P608*QL(K)))
230           T_PHY(I,K,J)=TL(K)
231           TH_PHY(I,K,J)=TL(K)*(1.E5/PLYR)**CAPA
232           P8W(I,K+1,J)=ETA1(K+1)*PDTOP+ETA2(K+1)*PDSL+PT
233           P_PHY(I,K,J)=PLYR
234           PI_PHY(I,K,J)=(PLYR*1.E-5)**CAPA
235           DZ(I,K,J)=TL(K)*(P608*QL(K)+1.)*R_D                           &
236      &                 *(P8W(I,K,J)-P8W(I,K+1,J))                       &
237      &                 /(P_PHY(I,K,J)*G)
238 !!!  &                 *ALOG(P8W(I,KFLIP,J)/P8W(I,KFLIP+1,J))/G         &
239 !!!  &                 *ALOG(PINT(I,K+1,J)/PINT(I,K,J))/G               &
240 !
241           RTHRATEN(I,K,J)=0.
242           THRATENLW(I,K,J)=0.
243           THRATENSW(I,K,J)=0.
244 !         PM2_5_DRY(I,K,J)=0.
245 !         PM2_5_WATER(I,K,J)=0.
246 
247         ENDDO
248 !
249         DO K=KTS+1,KTE
250           T8W(I,K,J)=0.5*(TL(K-1)+TL(K))
251         ENDDO
252         T8W(I,KTE+1,J)=-1.E20
253 !
254       ENDDO
255       ENDDO
256 !
257       ICLOUD=999
258 !
259       GMT=REAL(IHRST)
260 !
261 !-----------------------------------------------------------------------
262 !
263 !***  CALL THE INNER DRIVER.
264 !
265 !-----------------------------------------------------------------------
266 !
267       DO J=JMS,JME
268       DO K=KMS,KME
269       DO I=IMS,IME
270         QC1R(I,K,J)=0.
271         QI1R(I,K,J)=0.
272       ENDDO
273       ENDDO
274       ENDDO
275 !
276       DO J=MYJS2,MYJE2
277         DO K=KTS,KTE
278           DO I=MYIS1,MYIE1
279             QC1R(I,K,J)=MOIST(I,K,J,P_QC)
280             QI1R(I,K,J)=MOIST(I,K,J,P_QI)
281           ENDDO
282         ENDDO
283       ENDDO
284       DO J=JMS,JME
285         DO K=KMS,KME
286         DO I=IMS,IME
287           CLDFRA(I,K,J)=0.
288         ENDDO
289         ENDDO
290 !
291         DO I=IMS,IME
292           CFRACH(I,J)=0.
293           CFRACL(I,J)=0.
294           CFRACM(I,J)=0.
295           CZMEAN(I,J)=0.
296           SIGT4(I,J)=0.
297           TOTSWDN(I,J)=0.   ! TOTAL (clear+cloudy sky) shortwave down at the surface
298           TOTSWDNC(I,J)=0.  ! CLEAR SKY shortwave down at the surface
299           SWNETDN(I,J)=0.   ! Net (down - up) total (clear+cloudy sky) shortwave at the surface
300           TOTLWDN(I,J)=0.   ! Total longwave down at the surface
301           CUPPTR(I,J)=CUPPT(I,J)   ! Temporary array set to zero in radiation
302 !-- NOTE:  HBOTR, HTOPR are passed into radiation and set equal to HBOT, HTOP.  HBOT, HTOP are
303 !          reset to clear sky values to be used by the ARW.  At the bottom of this subroutine, 
304 !          HBOT, HTOP are re-defined again to values stored in HBOTR, HTOPR.  HBOT, HTOP are 
305 !          reset to clear sky values after the call to radiation and after the top of the hour
306 !          in subroutine CUCNVC below.
307         ENDDO
308       ENDDO
309 !
310       CALL SET_TILES(GRID,IDS+1,IDE-1,JDS+2,JDE-2,ITS,ITE,JTS,JTE)
311 !
312       CALL RADIATION_DRIVER(                                            &
313      &                  IDS=IDS,IDE=IDE,JDS=JDS,JDE=JDE,KDS=KDS,KDE=KDE &
314      &                 ,IMS=IMS,IME=IME,JMS=JMS,JME=JME,KMS=KMS,KME=KME &
315      &                 ,I_START=GRID%I_START,I_END=GRID%I_END           &
316      &                 ,J_START=GRID%J_START,J_END=GRID%J_END           &
317      &                 ,KTS=KTS,KTE=KTE,NUM_TILES=GRID%NUM_TILES        &
318      &                 ,ITIMESTEP=NTSD,DT=DT                            &
319      &                 ,cu_rad_feedback=config_flags%cu_rad_feedback    &
320 #ifdef WRF_CHEM
321      &        ,PM2_5_DRY=pm2_5_dry, PM2_5_WATER=pm2_5_water               &
322      &        ,PM2_5_DRY_EC=pm2_5_dry_ec                                  &
323      &        ,TAUAER300=tauaer1, TAUAER400=tauaer2, TAUAER600=tauaer3, TAUAER999=tauaer4 & ! jcb
324      &        ,GAER300=gaer1, GAER400=gaer2, GAER600=gaer3, GAER999=gaer4 & ! jcb
325      &        ,WAER300=waer1, WAER400=waer2, WAER600=waer3, WAER999=waer4 & ! jcb
326      &        ,qc_adjust=GD_CLOUD,qi_adjust=GD_CLOUD2                   &
327 #endif
328      &                 ,RTHRATENLW=THRATENLW,RTHRATENSW=THRATENSW       &
329      &                 ,RTHRATEN=RTHRATEN                               &
330      &                 ,GLW=TOTLWDN,GSW=SWNETDN,SWDOWN=TOTSWDN          &
331      &                 ,XLAT=XLAT,XLONG=XLON,ALBEDO=ALBEDO,EMISS=EPSR   &
332      &                 ,XICE=SICE,XLAND=XLAND,Z=Z,TSK=TSFC              &
333      &                 ,N_AEROSOLC=NUM_AEROSOLC,PAERLEV=GRID%PAERLEV    &
334      &                 ,CAM_ABS_DIM1=GRID%CAM_ABS_DIM1                  &
335      &                 ,CAM_ABS_DIM2=GRID%CAM_ABS_DIM2                  &
336      &                 ,CAM_ABS_FREQ_S=GRID%CAM_ABS_FREQ_S              &
337      &                 ,LEVSIZ=GRID%LEVSIZ,N_OZMIXM=NUM_OZMIXM          &
338      &                 ,HTOP=HTOP,HBOT=HBOT,CUPPT=CUPPTR                &
339      &                 ,HTOPR=HTOPR,HBOTR=HBOTR                         &
340      &                 ,VEGFRA=VEGFRC,SNOW=SNOW                         &
341      &                 ,RHO=RR,P8W=P8W,P=P_PHY,PI=PI_PHY                &
342      &                 ,DZ8W=DZ,T=T_PHY,T8W=T8W,GMT=GMT                 &
343      &                 ,JULDAY=JULDAY,JULYR=JULYR,NPHS=NPHS             &
344      &                 ,JULIAN=JULIAN,XTIME=XTIME                       &
345      &                 ,LW_PHYSICS=CONFIG_FLAGS%RA_LW_PHYSICS           &
346      &                 ,SW_PHYSICS=CONFIG_FLAGS%RA_SW_PHYSICS           &
347      &                 ,RADT=RADT,RA_CALL_OFFSET=GRID%RA_CALL_OFFSET    &
348      &                 ,STEPRA=NRAD,ICLOUD=ICLOUD                       &
349      &                 ,WARM_RAIN=WARM_RAIN                             & 
350      &                 ,SWDOWNC=TOTSWDNC,CLDFRA=CLFR                    &
351      &                 ,RSWTOA=RSWTOA,RLWTOA=RLWTOA                     &
352      &                 ,CZMEAN=CZMEAN,CFRACL=CFRACL                     &
353      &                 ,CFRACM=CFRACM,CFRACH=CFRACH                     &
354      &                 ,ACFRST=ACFRST,NCFRST=NCFRST                     &
355      &                 ,ACFRCV=ACFRCV,NCFRCV=NCFRCV                     &
356      &                 ,F_ICE_PHY=F_ICE,F_RAIN_PHY=F_RAIN               &
357      &                 ,QV=MOIST(IMS,KMS,JMS,P_QV),F_QV=F_QV            &
358      &                 ,QC=MOIST(IMS,KMS,JMS,P_QC),F_QC=F_QC                               &
359      &                 ,QR=MOIST(IMS,KMS,JMS,P_QR),F_QR=F_QR            &
360      &                 ,QI=MOIST(IMS,KMS,JMS,P_QI),F_QI=F_QI                               &
361      &                 ,QS=MOIST(IMS,KMS,JMS,P_QS),F_QS=F_QS            &
362      &                 ,QG=MOIST(IMS,KMS,JMS,P_QG),F_QG=F_QG          )
363 
364 !
365 !-----------------------------------------------------------------------
366 !
367 !***  UPDATE FLUXES AND TEMPERATURE TENDENCIES.
368 !
369 !-----------------------------------------------------------------------
370 !***  SHORTWAVE
371 !-----------------------------------------------------------------------
372 !
373 !-----------------------------------------------------------------------
374       IF(MOD(NTSD,NRADS)==0)THEN
375 !-----------------------------------------------------------------------
376 !
377         IF(CONFIG_FLAGS%RA_SW_PHYSICS/=GFDLSWSCHEME)THEN
378 !
379 !-----------------------------------------------------------------------
380 !***  COMPUTE CZMEAN FOR NON-GFDL SHORTWAVE
381 !-----------------------------------------------------------------------
382 !
383           DO J=MYJS,MYJE
384           DO I=MYIS,MYIE
385             CZMEAN(I,J)=0.
386             TOT(I,J)=0.
387           ENDDO
388           ENDDO
389 !
390           CALL CAL_MON_DAY(JULDAY,JULYR,JMONTH,JDAY)
391           IDAT(1)=JMONTH
392           IDAT(2)=JDAY
393           IDAT(3)=JULYR
394 !
395           DO II=0,NRADS,NPHS
396             TIMES=NTSD*DT+II*DT
397             CALL ZENITH(TIMES,DAYI,HOUR,IDAT,IHRST,GLON,GLAT,CZEN       &
398      &                 ,MYIS &
399      &                 ,MYIE &
400      &                 ,MYJS &
401      &                 ,MYJE &
402      &                 ,IDS,IDE,JDS,JDE,KDS,KDE                         &
403      &                 ,IMS,IME,JMS,JME,KMS,KME                         &
404      &                 ,ITS,ITE,JTS,JTE,KTS,KTE)
405             DO J=MYJS,MYJE
406             DO I=MYIS,MYIE
407               IF(CZEN(I,J)>0.)THEN
408                 CZMEAN(I,J)=CZMEAN(I,J)+CZEN(I,J)
409                 TOT(I,J)=TOT(I,J)+1.
410               ENDIF
411             ENDDO
412             ENDDO
413 !
414           ENDDO
415 !
416           DO J=MYJS,MYJE
417           DO I=MYIS,MYIE
418             IF(TOT(I,J)>0.)CZMEAN(I,J)=CZMEAN(I,J)/TOT(I,J)
419           ENDDO
420           ENDDO
421 !
422 !-----------------------------------------------------------------------
423 !***  COMPUTE TOTAL SFC SHORTWAVE DOWN FOR NON-GFDL SCHEMES
424 !-----------------------------------------------------------------------
425 !
426 !$omp parallel do                                                       &
427 !$omp& private(i,j)
428           DO J=MYJS2,MYJE2
429           DO I=MYIS1,MYIE1
430 !
431             IF(HBM2(I,J)>0.5)THEN
432               TOTSWDN(I,J)=SWNETDN(I,J)/(1.-ALBEDO(I,J))  
433 !--- No value currently available for clear-sky solar fluxes from
434 !    non GFDL schemes, though it's needed for air quality forecasts.
435 !    For the time being, set to the total downward solar fluxes.
436               TOTSWDNC(I,J)=TOTSWDN(I,J)
437             ENDIF
438 !
439           ENDDO
440           ENDDO
441 !
442         ENDIF   !End non-GFDL block
443 !-----------------------------------------------------------------------
444 !
445 !$omp parallel do                                                       &
446 !$omp& private(i,iendx,j,k)
447         DO J=MYJS2,MYJE2
448           IENDX=MYIE1
449           IF(MOD(J,2)==0.AND.ITE==IDE)IENDX=IENDX-1
450           DO I=MYIS1,IENDX
451 !
452             RSWIN(I,J)=TOTSWDN(I,J)
453             RSWINC(I,J)=TOTSWDNC(I,J)
454             RSWOUT(I,J)=TOTSWDN(I,J)-SWNETDN(I,J)
455 !
456             DO K=KTS,KTE
457               RSWTT(I,K,J)=THRATENSW(I,K,J)*PI_PHY(I,K,J)
458             ENDDO
459 !
460           ENDDO
461         ENDDO
462 !
463       ENDIF
464 !
465 !-----------------------------------------------------------------------
466 !***  LONGWAVE
467 !-----------------------------------------------------------------------
468 !
469       IF(MOD(NTSD,NRADL)==0)THEN
470 !
471 !$omp parallel do                                                       &
472 !$omp& private(i,iendx,j,k,lmhij)
473         DO J=MYJS2,MYJE2
474           IENDX=MYIE1
475           IF(MOD(J,2)==0.AND.ITE==IDE)IENDX=IENDX-1
476           DO I=MYIS1,IENDX
477 !
478             IF(HBM2(I,J)>0.5)THEN
479               LMHIJ=KTE+1-LMH(I,J)
480               TDUM=T(I,LMHIJ,J)
481               SIGT4(I,J)=STBOLT*TDUM*TDUM*TDUM*TDUM
482 !
483               DO K=KTS,KTE
484                 RLWTT(I,K,J)=THRATENLW(I,K,J)*PI_PHY(I,K,J)
485               ENDDO
486 !
487               RLWIN(I,J)=TOTLWDN(I,J)
488             ENDIF
489 !
490           ENDDO
491         ENDDO
492 !
493       ENDIF
494 !
495 !-- Store 3D cloud fractions & restore HBOT/HTOP arrays
496 !
497       DO J=MYJS2,MYJE2
498         IENDX=MYIE1
499         IF(MOD(J,2)==0.AND.ITE==IDE)IENDX=IENDX-1
500         DO K=KTS,KTE
501           DO I=MYIS1,IENDX
502             CLDFRA(I,K,J)=CLFR(I,K,J)
503           ENDDO
504         ENDDO
505         DO I=MYIS1,IENDX
506           HBOT(I,J)=HBOTR(I,J)
507           HTOP(I,J)=HTOPR(I,J)
508         ENDDO
509       ENDDO
510 !-----------------------------------------------------------------------
511 !***  ZERO OUT BOUNDARY ROWS.
512 !-----------------------------------------------------------------------
513 !
514       DO J=JTS,JTE
515       DO I=ITS,ITE
516         IF(HBM2(I,J)<0.5)THEN
517           ACFRST(I,J)=0.
518           ACFRCV(I,J)=0.
519           CFRACL(I,J)=0.
520           CFRACM(I,J)=0.
521           CFRACH(I,J)=0.
522           RSWTOA(I,J)=0.
523           RLWTOA(I,J)=0.
524         ENDIF
525       ENDDO
526       ENDDO
527 !
528 !-----------------------------------------------------------------------
529 !
530 
531       END SUBROUTINE RADIATION
532 !
533 !-----------------------------------------------------------------------
534 !***********************************************************************
535       SUBROUTINE TURBL(NTSD,DT,NPHS,RESTRT                              &
536      &                ,N_MOIST,NSOIL,SLDPTH,DZSOIL                      &
537      &                ,DETA1,DETA2,AETA1,AETA2,ETA1,ETA2,PDTOP,PT       &
538      &                ,SM,LMH,HTM,VTM,HBM2,VBM2,DX_ARRAY,DFRLG          &
539      &                ,CZEN,CZMEAN,SIGT4,RLWIN,RSWIN,RADOT              &
540 !- RLWIN/RSWIN - downward longwave/shortwave at the surface (also TOTLWDN/TOTSWDN in RADIATION)
541      &                ,PD,RES,PINT,T,Q,CWM,F_ICE,F_RAIN,SR              &
542      &                ,Q2,U,V,THS,TSFC,SST,PREC,SNO,ZERO_3D             &
543      &                ,FIS,Z0,Z0BASE,USTAR,PBLH,LPBL,EL_MYJ             &
544      &                ,MOIST,RMOL                                       &
545      &                ,EXCH_H,AKHS,AKMS,AKHS_OUT,AKMS_OUT               &
546      &                ,THZ0,QZ0,UZ0,VZ0,QS,MAVAIL                       &
547      &                ,STC,SMC,CMC,SMSTAV,SMSTOT,SSROFF,BGROFF          &
548      &                ,IVGTYP,ISLTYP,VEGFRC,SHDMIN,SHDMAX,GRNFLX        &
549      &                ,SFCEXC,ACSNOW,ACSNOM,SNOPCX,SICE,TG,SOILTB       &
550      &                ,ALBASE,MXSNAL,ALBEDO,SH2O,SI,EPSR                &
551      &                ,U10,V10,TH10,Q10,TSHLTR,QSHLTR,PSHLTR            &
552      &                ,T2,QSG,QVG,QCG,SOILT1,TSNAV,SMFR3D,KEEPFR3DFLAG  &
553      &                ,TWBS,QWBS,SFCSHX,SFCLHX,SFCEVP                   &
554      &                ,POTEVP,POTFLX,SUBSHX                             &
555      &                ,APHTIM,ARDSW,ARDLW,ASRFC                         &
556      &                ,RSWOUT,RSWTOA,RLWTOA                             &
557      &                ,ASWIN,ASWOUT,ASWTOA,ALWIN,ALWOUT,ALWTOA          &
558      &                ,UZ0H,VZ0H,DUDT,DVDT                              & 
559      &                ,RTHBLTEN,RQVBLTEN                                & 
560      &                ,PCPFLG,DDATA                                     & ! PRECIP ASSIM
561      &                ,GRID,CONFIG_FLAGS                                &
562      &                ,IHE,IHW,IVE,IVW                                  &
563      &                ,IDS,IDE,JDS,JDE,KDS,KDE                          &
564      &                ,IMS,IME,JMS,JME,KMS,KME                          &
565      &                ,ITS,ITE,JTS,JTE,KTS,KTE)
566 !***********************************************************************
567 !$$$  SUBPROGRAM DOCUMENTATION BLOCK
568 !                .      .    .     
569 ! SUBPROGRAM:    TURBL       TURBULENCE OUTER DRIVER
570 !   PRGRMMR: BLACK           ORG: W/NP22     DATE: 02-04-19       
571 !     
572 ! ABSTRACT:
573 !     TURBL DRIVES THE TURBULENCE SCHEMES
574 !     
575 ! PROGRAM HISTORY LOG (with changes to called routines) :
576 !   95-03-15  JANJIC     - ORIGINATOR OF THE SUBROUTINES CALLED
577 !   BLACK & JANJIC       - ORIGINATORS OF THE DRIVER
578 !   95-03-28  BLACK      - CONVERSION FROM 1-D TO 2-D IN HORIZONTAL
579 !   96-03-29  BLACK      - ADDED EXTERNAL EDGE; REMOVED SCRCH COMMON
580 !   96-07-19  MESINGER   - ADDED Z0 EFFECTIVE
581 !   98-??-??  TUCCILLO   - MODIFIED FOR CLASS VIII PARALLELISM
582 !   98-10-27  BLACK      - PARALLEL CHANGES INTO MOST RECENT CODE
583 !   02-01-10  JANJIC     - MOIST TURBULENCE (DRIVER, MIXLEN, VDIFH)
584 !   02-01-10  JANJIC     - VERT. DIF OF Q2 INCREASED (Grenier & Bretherton)
585 !   02-02-02  JANJIC     - NEW SFCDIF
586 !   02-04-19  BLACK      - ORIGINATOR OF THIS OUTER DRIVER FOR WRF
587 !   02-05-03  JANJIC     - REMOVAL OF SUPERSATURATION AT 2m AND 10m
588 !   04-11-18  BLACK      - THREADED
589 !     
590 ! USAGE: CALL TURBL FROM SOLVE_NMM
591 !
592 ! ATTRIBUTES:
593 !   LANGUAGE: FORTRAN 90
594 !   MACHINE : IBM
595 !$$$  
596 !-----------------------------------------------------------------------
597 !
598       IMPLICIT NONE
599 !
600 !-----------------------------------------------------------------------
601 !
602       INTEGER,INTENT(IN) :: IDS,IDE,JDS,JDE,KDS,KDE                     &
603      &                     ,IMS,IME,JMS,JME,KMS,KME                     &
604      &                     ,ITS,ITE,JTS,JTE,KTS,KTE                     &
605      &                     ,N_MOIST,NPHS,NSOIL,NTSD
606 !
607       INTEGER, DIMENSION(JMS:JME),INTENT(IN) :: IHE,IHW,IVE,IVW
608 !
609       INTEGER,DIMENSION(IMS:IME,JMS:JME),INTENT(IN) :: ISLTYP,IVGTYP    &
610      &                                                ,LMH
611 !
612       INTEGER,DIMENSION(IMS:IME,JMS:JME),INTENT(OUT) :: LPBL
613 !
614       REAL,INTENT(IN) :: DT,PDTOP,PT
615 !
616       REAL,INTENT(INOUT) :: APHTIM,ARDSW,ARDLW,ASRFC
617 !
618       REAL,DIMENSION(KMS:KME-1),INTENT(IN) :: AETA1,AETA2,DETA1,DETA2
619 !
620       REAL,DIMENSION(KMS:KME),INTENT(IN) :: DFRLG,ETA1,ETA2
621 !
622       REAL,DIMENSION(IMS:IME,JMS:JME),INTENT(INOUT) :: ALBASE,MXSNAL
623 !
624       REAL,DIMENSION(IMS:IME,JMS:JME),INTENT(IN) :: CZEN,CZMEAN         &
625      &                                             ,DX_ARRAY            &
626      &                                             ,FIS,HBM2            &
627      &                                             ,PD,RES              &
628      &                                             ,RLWIN,RLWTOA        &
629      &                                             ,RSWIN,RSWOUT,RSWTOA &
630      &                                             ,SHDMIN,SHDMAX       &
631 !    &                                             ,SICE,SIGT4,SM,SR    & !Bandaid
632      &                                             ,SICE,SIGT4          &
633      &                                             ,SST,TG,VBM2,VEGFRC
634 !
635       REAL,DIMENSION(IMS:IME,JMS:JME),INTENT(INOUT) :: SM,EPSR,SR         !Bandaid
636 !
637       REAL,DIMENSION(IMS:IME,JMS:JME),INTENT(OUT) :: GRNFLX,QWBS,RADOT  &
638                                                     ,SFCEXC,SMSTAV      &
639                                                     ,SOILTB,TWBS
640 !
641       REAL,DIMENSION(IMS:IME,JMS:JME),INTENT(INOUT) :: ACSNOM,ACSNOW    &
642      &                                                ,AKHS,AKMS        &
643      &                                                ,ALBEDO           &
644      &                                                ,MAVAIL           &
645      &                                                ,BGROFF,CMC       &
646      &                                                ,PBLH,POTEVP      &
647      &                                                ,POTFLX,PREC      &
648      &                                                ,QCG,QS,QSG       &
649      &                                                ,QVG,QZ0          &
650      &                                                ,SFCEVP           &
651      &                                                ,SFCLHX,SFCSHX    &
652      &                                                ,SI,SMSTOT        &
653      &                                                ,SNO,SNOPCX       &
654      &                                                ,SOILT1           &
655      &                                                ,SSROFF,SUBSHX    &
656      &                                                ,T2,THS,THZ0      &
657      &                                                ,TSFC,TSNAV       &
658      &                                                ,USTAR,UZ0,UZ0H   &
659      &                                                ,VZ0,VZ0H         &
660      &                                                ,Z0,Z0BASE
661 !
662       REAL,DIMENSION(IMS:IME,JMS:JME),INTENT(OUT) :: AKHS_OUT,AKMS_OUT  &
663      &                                              ,ALWIN,ALWOUT       &
664      &                                              ,ALWTOA,ASWIN       &
665      &                                              ,ASWOUT,ASWTOA      &
666      &                                              ,PSHLTR,Q10,QSHLTR  &
667      &                                              ,TH10,TSHLTR        &
668      &                                              ,U10,V10
669 !
670       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(INOUT) :: CWM      &
671      &                                                        ,DUDT     &
672      &                                                        ,DVDT     &
673      &                                                        ,EXCH_H   &
674      &                                                        ,F_ICE    &
675      &                                                        ,F_RAIN   &
676      &                                                        ,Q,Q2     &
677      &                                                        ,T,U,V
678       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(INOUT) :: RQVBLTEN,RTHBLTEn
679       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME,N_MOIST),INTENT(INOUT) :: MOIST
680 !
681       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(IN) :: HTM,VTM
682       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(IN) :: PINT
683 !
684       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(INOUT) :: ZERO_3D
685 !
686       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(OUT) :: EL_MYJ
687 !
688       REAL,DIMENSION(NSOIL),INTENT(IN) :: DZSOIL,SLDPTH
689 !
690       REAL,DIMENSION(IMS:IME,NSOIL,JMS:JME),INTENT(INOUT) :: KEEPFR3DFLAG &
691      &                                                      ,SH2O,SMC     &
692      &                                                      ,SMFR3D,STC
693 !
694       LOGICAL,INTENT(IN) :: RESTRT
695 !
696       TYPE(DOMAIN),TARGET :: GRID
697 !
698       TYPE(GRID_CONFIG_REC_TYPE),INTENT(IN) :: CONFIG_FLAGS
699 !
700 !  For precip assimilation:
701       LOGICAL,INTENT(IN) :: PCPFLG
702       REAL,DIMENSION(IMS:IME,JMS:JME),INTENT(IN) :: DDATA
703 !
704 !-----------------------------------------------------------------------
705 !***
706 !***  LOCAL VARIABLES
707 !***
708 !-----------------------------------------------------------------------
709       INTEGER :: I,I_M,IDUMMY,IEND,ISFFLX,ISTR,J,K,KOUNT_ALL,LENGTH_ROW &
710      &          ,LLIJ,LLMH,LLYR,N,SST_UPDATE
711 !
712       INTEGER,DIMENSION(IMS:IME,JMS:JME) :: KPBL,LOWLYR
713 !
714       REAL :: TRESH=0.95
715 !
716       REAL :: ALTITUDE,CWML,DQDT,DTDT,DTPHS,DX,DZHALF,FACTR,FACTRL      &
717      &       ,G_INV,PDSL,PLYR,PSFC,QI,QL,QOLD,QR,QW,RATIOMX,RDTPHS      &
718      &       ,ROG,RWMSK,SDEPTH,SNO_FACTR,TL,TLMH,TLMH4,TNEW,TSFC2       &
719      &       ,U_FRAME,V_FRAME,WMSK,XLVRW
720 !
721       REAL :: APES,CKLQ,FACTOR,FFS,PQ0X,Q2SAT,QFC1,QLOWX,RLIVWV         &
722      &       ,THBOT
723 !
724       REAL,DIMENSION(IMS:IME,JMS:JME) :: BR,CHKLOWQ,CT,CWMLOW,ELFLX     &
725      &                                  ,EXNSFC,FACTRS,FLHC,FLQC,GZ1OZ0 &
726      &                                  ,ONE,PLM,PSFC_OUT,PSIH,PSIM     &
727      &                                  ,Q2X,QLOW,RAIN,RAINBL           &
728      &                                  ,RLW_DN_SFC,RMOL,RSW_NET_SFC    &
729      &                                  ,RSW_DN_SFC                     &
730      &                                  ,SFCEVPX,SFCZ,SNOW,SNOWC,SNOWH  &
731      &                                  ,TH2X,THLOW,TLOW,VGFRCK         &
732      &                                  ,WSPD,XLAND,ZERO_2D,EMISS
733 !
734       REAL,DIMENSION(IMS:IME,KMS:KME-1,JMS:JME) :: EXNER
735 !
736       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME) :: DZ,P8W                 &
737      &                                          ,P_PHY,PI_PHY           &
738      &                                          ,RQCBLTEN,RQIBLTEN      &
739      &                                          ,RR   &
740 !    &                                          ,RQVBLTEN,RR,RTHBLTEN   &
741      &                                          ,T_PHY,TH_PHY,TKE       &
742      &                                          ,U_PHY,V_PHY,Z
743 !
744       REAL,DIMENSION(IMS:IME,NSOIL,JMS:JME) :: ZERO_SOIL
745 !
746       LOGICAL :: E_BDY,WARM_RAIN
747 !
748       INTEGER :: ucmcall
749 !
750 !-----------------------------------------------------------------------
751 !-----------------------------------------------------------------------
752      ucmcall=config_flags%ucmcall
753 !
754       DTPHS=NPHS*DT
755       RDTPHS=1./DTPHS
756       G_INV=1./G
757       ROG=R_D*G_INV
758       FACTOR=-XLV*RHOWATER/DTPHS
759 !
760       U_FRAME=0.
761       V_FRAME=0.
762 !
763       IDUMMY=0
764       ISFFLX=1
765       DX=0.
766       SST_UPDATE=0
767 !
768       DO J=JMS,JME
769       DO I=IMS,IME
770         UZ0H(I,J)=0.
771         VZ0H(I,J)=0.
772         ONE(I,J)=1.
773         RMOL(I,J)=0.  !Reciprocal of Monin-Obukhov length
774         SFCEVPX(I,J)=0.  !Dummy for accumulated latent energy, not flux
775       ENDDO
776       ENDDO
777 !
778       IF(MODEL_CONFIG_REC%SF_SURFACE_PHYSICS(GRID%ID)==99)THEN
779         SNO_FACTR=1.
780       ELSE
781         SNO_FACTR=0.001
782       ENDIF
783 !
784 !$omp parallel do                                                       &
785 !$omp& private(i,j)
786       DO J=MYJS,MYJE
787       DO I=MYIS,MYIE
788         LOWLYR(I,J)=1
789         VGFRCK(I,J)=100.*VEGFRC(I,J)
790         SNOW(I,J)=SNO(I,J)
791         SNOWH(I,J)=SI(I,J)*SNO_FACTR
792         XLAND(I,J)=SM(I,J)+1.
793         T2(I,J)=TSFC(I,J)
794         EMISS(I,J)=EPSR(I,J)
795       ENDDO
796       ENDDO
797 !
798       IF(NTSD==0)THEN
799 !$omp parallel do                                                       &
800 !$omp& private(i,j)
801         DO J=MYJS,MYJE
802         DO I=MYIS,MYIE
803           Z0BASE(I,J)=Z0(I,J)
804           IF(SM(I,J)>0.5.AND.SICE(I,J)>0.5)THEN  !Bandaid
805             SM(I,J)=0.        
806           ENDIF              
807         ENDDO
808         ENDDO
809       ENDIF
810 !
811 !$omp parallel do                                                       &
812 !$omp& private(i,j,k)
813       DO J=MYJS,MYJE
814       DO K=KTS,KTE+1
815       DO I=MYIS,MYIE
816         Z(I,K,J)=0.
817         DZ(I,K,J)=0.
818         EXCH_H(I,K,J)=0.
819       ENDDO
820       ENDDO
821       ENDDO
822 !
823 !-----------------------------------------------------------------------
824 !
825 !***  PREPARE NEEDED ARRAYS
826 !
827 !-----------------------------------------------------------------------
828 !
829 !$omp parallel do                                                       &
830 !$omp& private(cwml,factrl,i,j,k,llij,llmh,pdsl,plyr,psfc,qi,ql,qr,qw   &
831 !$omp&        ,tl,tlmh,tlmh4)
832       DO J=MYJS,MYJE
833       DO I=MYIS,MYIE
834 !
835         LLMH=LMH(I,J)
836         PDSL=PD(I,J)*RES(I,J)
837 !!!     PSFC=PD(I,J)+PDTOP+PT
838 !!!     P8W(I,KTS,J)=PSFC
839         P8W(I,KTS,J)=PINT(I,KTS,J)
840         PSFC=PINT(I,KTS,J)
841         LOWLYR(I,J)=KTE+1-LLMH
842         EXNSFC(I,J)=(1.E5/PSFC)**CAPA
843         THS(I,J)=(SST(I,J)*EXNSFC(I,J))*SM(I,J)+THS(I,J)*(1.-SM(I,J))
844         TSFC(I,J)=THS(I,J)/EXNSFC(I,J)
845         SFCZ(I,J)=FIS(I,J)*G_INV
846         ZERO_2D(I,J)=0.
847 !YL     RAIN(I,J)=PREC(I,J)*RHOWATER
848         IF (PCPFLG.AND.DDATA(I,J)<100.)THEN
849           RAIN(I,J)=DDATA(I,J)*RHOWATER
850         ELSE
851           RAIN(I,J)=PREC(I,J)*RHOWATER
852         ENDIF
853 !YL
854         RAINBL(I,J)=0.
855         IF(SNO(I,J)>0.)SNOWC(I,J)=1.
856         LLIJ=LOWLYR(I,J)
857         PLM(I,J)=(PINT(I,LLIJ,J)+PINT(I,LLIJ+1,J))*0.5
858         TH2X(I,J)=T(I,LLIJ,J)*(1.E5/PLM(I,J))**CAPA
859         Q2X(I,J)=Q(I,LLIJ,J)
860 !
861 !-----------------------------------------------------------------------
862 !*** LONG AND SHORTWAVE FLUX AT GROUND SURFACE
863 !-----------------------------------------------------------------------
864 !
865         IF(CZMEAN(I,J)>0.)THEN
866           FACTRS(I,J)=CZEN(I,J)/CZMEAN(I,J)
867         ELSE
868           FACTRS(I,J)=0.
869         ENDIF
870 !
871         IF(SIGT4(I,J)>0.)THEN
872           TLMH=T(I,LLIJ,J)
873           FACTRL=STBOLT*TLMH*TLMH*TLMH*TLMH/SIGT4(I,J)
874         ELSE
875           FACTRL=0.
876         ENDIF
877 !     
878 !- RLWIN/RSWIN - downward longwave/shortwave at the surface
879 !
880         RLW_DN_SFC(I,J)=RLWIN(I,J)*HBM2(I,J)*FACTRL
881         RSW_NET_SFC(I,J)=(RSWIN(I,J)-RSWOUT(I,J))*HBM2(I,J)*FACTRS(I,J)
882 !
883 !- Instant downward solar for nmm_lsm
884 !
885         RSW_DN_SFC(I,J)=RSWIN(I,J)*HBM2(I,J)*FACTRS(I,J)
886 !
887 !-----------------------------------------------------------------------
888 !***  FILL THE ARRAYS FOR CALLING THE INNER DRIVER.
889 !-----------------------------------------------------------------------
890 !
891         Z(I,KTS,J)=SFCZ(I,J)
892 !
893         DO K=KTS,KTE
894           Q2(I,K,J)=AMAX1(Q2(I,K,J)*HBM2(I,J),EPSQ2)
895           QL=AMAX1(Q(I,K,J),EPSQ)
896           PLYR=(PINT(I,K,J)+PINT(I,K+1,J))*0.5
897 !!!       PLYR=AETA1(K)*PDTOP+AETA2(K)*PDSL+PT
898           TL=T(I,K,J)
899           CWML=CWM(I,K,J)
900 !
901           RR(I,K,J)=PLYR/(R_D*TL)
902           T_PHY(I,K,J)=TL
903 !
904           EXNER(I,K,J)=(1.E5/PLYR)**CAPA
905           PI_PHY(I,K,J)=1./EXNER(I,K,J)
906           TH_PHY(I,K,J)=TL*EXNER(I,K,J)
907           P8W(I,K+1,J)=PINT(I,K+1,J)
908 !!!       P8W(I,K+1,J)=ETA1(K+1)*PDTOP+ETA2(K+1)*PDSL+PT
909           P_PHY(I,K,J)=PLYR
910           TKE(I,K,J)=0.5*Q2(I,K,J)
911 !
912           RTHBLTEN(I,K,J)=0.
913           RQVBLTEN(I,K,J)=0.
914           RQCBLTEN(I,K,J)=0.
915           RQIBLTEN(I,K,J)=0.
916 !
917           Z(I,K+1,J)=Z(I,K,J)+TL/PLYR                                   &
918      &                  *(DETA1(K)*PDTOP+DETA2(K)*PDSL)*ROG             &
919                         *(Q(I,K,J)*P608-CWML+1.)
920           Z(I,K+1,J)=(Z(I,K+1,J)-DFRLG(K+1))*HTM(I,K,J)+DFRLG(K+1)
921 !!!       FACTR=1.-HTM(I,K,J)
922 !!!       Z(I,K+1,J)=Z(I,K+1,J)*HTM(I,K,J)+FACTR*DFRLG(K+1)
923           DZ(I,K,J)=Z(I,K+1,J)-Z(I,K,J)
924         ENDDO
925       ENDDO
926       ENDDO
927 !
928 !$omp parallel do                                                       &
929 !$omp& private(i,j,llyr,qlowx)
930       DO J=MYJS,MYJE
931       DO I=MYIS,MYIE
932         TWBS(I,J)=0.
933         QWBS(I,J)=0.
934         LLYR=LOWLYR(I,J)
935         THLOW(I,J)=TH_PHY(I,LLYR,J)
936         TLOW(I,J)=T_PHY(I,LLYR,J)
937         QLOW(I,J)=MAX(Q(I,LLYR,J),EPSQ)
938         QLOWX=QLOW(I,J)/(1.-QLOW(I,J))
939         QLOW(I,J)=QLOWX/(1.+QLOWX)
940         CWMLOW(I,J)=CWM(I,LLYR,J)
941         PBLH(I,J)=MAX(PBLH(I,J),0.)
942         PBLH(I,J)=MIN(PBLH(I,J),Z(I,KTE,J))
943       ENDDO
944       ENDDO
945 !-----------------------------------------------------------------------
946 !
947 !***  COMPUTE VELOCITY COMPONENTS AT MASS POINTS
948 !
949 !-----------------------------------------------------------------------
950 !$omp parallel do                                                       &
951 !$omp& private(i,j,k,rwmsk,wmsk)
952       DO J=MYJS1_P1,MYJE1_P1
953 !
954         DO K=KTS,KTE
955           DO I=MYIS_P1,MYIE_P1
956             WMSK=VTM(I+IHE(J),K,J)+VTM(I+IHW(J),K,J)                    &
957      &          +VTM(I,K,J+1)+VTM(I,K,J-1)
958             IF(WMSK>0.)THEN
959               RWMSK=1./WMSK
960               U_PHY(I,K,J)=(U(I+IHE(J),K,J)*VTM(I+IHE(J),K,J)           &
961      &                         +U(I+IHW(J),K,J)*VTM(I+IHW(J),K,J)       &
962      &                         +U(I,K,J+1)*VTM(I,K,J+1)                 &
963      &                         +U(I,K,J-1)*VTM(I,K,J-1))*RWMSK
964               V_PHY(I,K,J)=(V(I+IHE(J),K,J)*VTM(I+IHE(J),K,J)           &
965      &                         +V(I+IHW(J),K,J)*VTM(I+IHW(J),K,J)       &
966      &                         +V(I,K,J+1)*VTM(I,K,J+1)                 &
967      &                         +V(I,K,J-1)*VTM(I,K,J-1))*RWMSK
968             ELSE
969               U_PHY(I,K,J)=0.
970               V_PHY(I,K,J)=0.
971             ENDIF
972           ENDDO
973         ENDDO
974       ENDDO
975 !
976 !$omp parallel do                                                       &
977 !$omp& private(i,iend,istr,j)
978       DO J=MYJS1_P1,MYJE1_P1
979         IF(MOD(J,2)==0)THEN
980           ISTR=MYIS_P1
981           IEND=MIN(MYIE_P1,IDE-1)
982         ELSE
983           ISTR=MAX(MYIS_P1,IDS+1)
984           IEND=MIN(MYIE_P1,IDE-1)
985         ENDIF
986 !     
987         DO I=ISTR,IEND
988           UZ0H(I,J)=(UZ0(I+IHE(J),J)+UZ0(I+IHW(J),J)                    &
989      &              +UZ0(I,J+1)+UZ0(I,J-1))*0.25
990 !!!  &              +UZ0(I,J+1)+UZ0(I,J-1))*HBM2(I,J)*0.25
991           VZ0H(I,J)=(VZ0(I+IHE(J),J)+VZ0(I+IHW(J),J)                    &
992      &              +VZ0(I,J+1)+VZ0(I,J-1))*0.25
993 !!!  &              +VZ0(I,J+1)+VZ0(I,J-1))*HBM2(I,J)*0.25
994         ENDDO
995       ENDDO
996 !-----------------------------------------------------------------------
997 !
998 !***  CALL SURFACE LAYER AND LAND SURFACE PHYSICS
999 !
1000 !-----------------------------------------------------------------------
1001 !
1002       CALL SET_TILES(GRID,IDS,IDE-1,JDS+1,JDE-1,ITS,ITE,JTS,JTE)
1003 !
1004       DO J=JTS,JTE  !jm was JTS
1005       DO I=ITS,ITE
1006         IF(MODEL_CONFIG_REC%SF_SURFACE_PHYSICS(GRID%ID)==99)THEN
1007           ONE(I,J)=1.
1008         ELSE
1009 !tgs  -  MAVAIL should not be equal to 1. for other LSMs
1010           ONE(I,J)=MAVAIL(I,J)
1011         ENDIF 
1012       ENDDO
1013       ENDDO
1014 !
1015       CALL SURFACE_DRIVER(                                              &
1016      &           ACSNOM=ACSNOM,ACSNOW=ACSNOW,AKHS=AKHS,AKMS=AKMS        &
1017      &          ,ALBEDO=ALBEDO,BR=BR,CANWAT=CMC,CHKLOWQ=CHKLOWQ         &
1018      &          ,DT=DT,DX=DX,DZ8W=DZ,DZS=DZSOIL,GLW=RLW_DN_SFC          &
1019      &          ,GRDFLX=GRNFLX,GSW=RSW_NET_SFC,SWDOWN=RSW_DN_SFC        &
1020      &          ,GZ1OZ0=GZ1OZ0,HFX=TWBS                                 &
1021      &          ,HT=SFCZ,IFSNOW=IDUMMY,ISFFLX=ISFFLX,ISLTYP=ISLTYP      &
1022      &          ,ITIMESTEP=NTSD,IVGTYP=IVGTYP,LOWLYR=LOWLYR             &
1023      &          ,MAVAIL=ONE,RMOL=RMOL,NUM_SOIL_LAYERS=NSOIL,P8W=P8W &
1024      &          ,PBLH=PBLH,PI_PHY=PI_PHY,PSHLTR=PSHLTR,PSIH=PSIH        &
1025      &          ,PSIM=PSIM,P_PHY=P_PHY,Q10=Q10,Q2=Q2X,QFX=QWBS,QSFC=QS  &
1026      &          ,QSHLTR=QSHLTR,QZ0=QZ0,RAINCV=RAIN                      &
1027      &          ,RHO=RR,SFCEVP=SFCEVPX,SFCEXC=SFCEXC,SFCRUNOFF=SSROFF   &
1028      &          ,SMOIS=SMC,SMSTAV=SMSTAV,SMSTOT=SMSTOT,SNOALB=MXSNAL    &
1029      &          ,SNOW=SNOW,SNOWC=SNOWC,SNOWH=SNOWH,STEPBL=NPHS          &
1030      &          ,SST=SST,SST_UPDATE=SST_UPDATE                          &
1031      &          ,TH10=TH10,TH2=TH2X,T2=T2,THZ0=THZ0,TH_PHY=TH_PHY       &
1032      &          ,TMN=TG,TSHLTR=TSHLTR,TSK=TSFC,TSLB=STC,T_PHY=T_PHY     &
1033      &          ,U10=U10,UDRUNOFF=BGROFF,UST=USTAR,UZ0=UZ0H             &
1034      &          ,U_FRAME=U_FRAME,U_PHY=U_PHY,V10=V10,VEGFRA=VGFRCK      &
1035      &          ,VZ0=VZ0H,V_FRAME=V_FRAME,V_PHY=V_PHY                   &
1036      &          ,WARM_RAIN=WARM_RAIN,WSPD=WSPD,XICE=SICE                &
1037      &          ,XLAND=XLAND,Z=Z,ZNT=Z0,ZS=SLDPTH,CT=CT,TKE_MYJ=TKE     &
1038      &          ,ALBBCK=ALBASE,LH=ELFLX,SH2O=SH2O,SHDMAX=SHDMAX         &
1039      &          ,SHDMIN=SHDMIN,Z0=Z0BASE,FLQC=FLQC,FLHC=FLHC            &
1040      &          ,PSFC=PSFC_OUT,EMISS=EPSR                               &
1041      &          ,SF_SFCLAY_PHYSICS=CONFIG_FLAGS%SF_SFCLAY_PHYSICS       &
1042      &          ,SF_SURFACE_PHYSICS=CONFIG_FLAGS%SF_SURFACE_PHYSICS     &
1043      &          ,RA_LW_PHYSICS=CONFIG_FLAGS%RA_LW_PHYSICS               &
1044      &          ,UCMCALL=ucmcall                                        &
1045      &          ,IDS=IDS,IDE=IDE,JDS=JDS,JDE=JDE,KDS=KDS,KDE=KDE        &
1046      &          ,IMS=IMS,IME=IME,JMS=JMS,JME=JME,KMS=KMS,KME=KME        &
1047      &          ,I_START=GRID%I_START,I_END=GRID%I_END                  &
1048      &          ,J_START=GRID%J_START,J_END=GRID%J_END                  &
1049      &          ,KTS=KTS,KTE=KTE,NUM_TILES=GRID%NUM_TILES               &
1050            ! Optional args
1051      &          ,QV_CURR=MOIST(IMS,KMS,JMS,P_QV),F_QV=F_QV              &
1052      &          ,QC_CURR=MOIST(IMS,KMS,JMS,P_QC),F_QC=F_QC              &
1053      &          ,QR_CURR=MOIST(IMS,KMS,JMS,P_QR),F_QR=F_QR              &
1054      &          ,QI_CURR=MOIST(IMS,KMS,JMS,P_QI),F_QI=F_QI              &
1055      &          ,QS_CURR=MOIST(IMS,KMS,JMS,P_QS),F_QS=F_QS              & 
1056      &          ,QG_CURR=MOIST(IMS,KMS,JMS,P_QG),F_QG=F_QG              &
1057      &          ,RAINBL=RAINBL                                          &
1058 ! for RUCLSM
1059      &          ,QSG=QSG, QVG=QVG, QCG=QCG, SOILT1=SOILT1               &
1060      &          ,TSNAV=TSNAV, SMFR3D=SMFR3D, KEEPFR3DFLAG=KEEPFR3DFLAG  &
1061      &          ,POTEVP=POTEVP,SNOPCX=SNOPCX,SOILTB=SOILTB,SR=SR)
1062 !
1063 !-----------------------------------------------------------------------
1064 !
1065 !***  CALL FREE ATMOSPHERE TURBULENCE
1066 !
1067 !-----------------------------------------------------------------------
1068 !
1069 !$omp parallel do                                                       &
1070 !$omp& private(i,j,k)
1071       DO J=JMS,JME
1072       DO K=KMS,KME
1073       DO I=IMS,IME
1074         DUDT(I,K,J)=0.
1075         DVDT(I,K,J)=0.
1076       ENDDO
1077       ENDDO
1078       ENDDO
1079 !
1080 !***  THE SURFACE EXCHANGE COEFFICIENTS AKHS AND AKMS ARE ACTUALLY
1081 !***  MULTIPLIED BY HALF THE DEPTH OF THE LOWEST LAYER.  WE MUST RETAIN
1082 !***  THOSE VALUES FOR THE NEXT TIMESTEP SO USE AUXILLIARY ARRAYS FOR
1083 !***  THE OUTPUT.
1084 !
1085 !$omp parallel do                                                       &
1086 !$omp& private(dzhalf,i,j)
1087       DO J=JTS,JTE
1088       DO I=ITS,ITE
1089         DZHALF=0.5*DZ(I,KTS,J)
1090         AKHS_OUT(I,J)=AKHS(I,J)*DZHALF
1091         AKMS_OUT(I,J)=AKMS(I,J)*DZHALF
1092       ENDDO
1093       ENDDO
1094 !
1095       CALL PBL_DRIVER(                                                &
1096      &                ITIMESTEP=NTSD,DT=DT                            &
1097      &               ,U_FRAME=U_FRAME,V_FRAME=V_FRAME                 &
1098      &               ,RUBLTEN=DUDT,RVBLTEN=DVDT,RTHBLTEN=RTHBLTEN     &
1099      &               ,RQVBLTEN=RQVBLTEN,RQCBLTEN=RQCBLTEN             &
1100      &               ,RQSBLTEN=RQIBLTEN                               &
1101      &               ,TSK=TSFC,XLAND=XLAND,ZNT=Z0,HT=SFCZ             &
1102      &               ,UST=USTAR, PBLH=PBLH                            &
1103      &               ,HFX=TWBS,QFX=QWBS,  GRDFLX=GRNFLX               &
1104      &               ,U_PHY=U_PHY,V_PHY=V_PHY,TH_PHY=TH_PHY,RHO=RR    &
1105      &               ,P_PHY=P_PHY,PI_PHY=PI_PHY,P8W=P8W,T_PHY=T_PHY   &
1106      &               ,DZ8W=DZ,Z=Z,TKE_MYJ=TKE,EL_MYJ=EL_MYJ           &
1107      &               ,EXCH_H=EXCH_H,AKHS=AKHS,AKMS=AKMS               &
1108      &               ,THZ0=THZ0,QZ0=QZ0,UZ0=UZ0H,VZ0=VZ0H             &
1109      &               ,QSFC=QS,LOWLYR=LOWLYR                           &
1110      &               ,PSIM=PSIM,PSIH=PSIH,GZ1OZ0=GZ1OZ0               &
1111      &               ,WSPD=WSPD,BR=BR,CHKLOWQ=CHKLOWQ                 &
1112      &               ,DX=DX,STEPBL=NPHS,WARM_RAIN=WARM_RAIN           &
1113      &               ,KPBL=KPBL,CT=CT,LH=ELFLX,SNOW=SNOW,XICE=SICE    &
1114      &               ,BL_PBL_PHYSICS=config_flags%bl_pbl_physics      &
1115      &               ,RA_LW_PHYSICS=config_flags%ra_lw_physics        &
1116      &               ,IDS=IDS,IDE=IDE,JDS=JDS,JDE=JDE,KDS=KDS,KDE=KDE &
1117      &               ,IMS=IMS,IME=IME,JMS=JMS,JME=JME,KMS=KMS,KME=KME &
1118      &               ,I_START=GRID%I_START,I_END=GRID%I_END           &
1119      &               ,J_START=GRID%J_START,J_END=GRID%J_END           &
1120      &               ,KTS=KTS,KTE=KTE,NUM_TILES=GRID%NUM_TILES        &
1121                 ! Optional args
1122      &               ,QV_CURR=moist(IMS,KMS,JMS,P_QV) , F_QV=F_QV     &
1123      &               ,QC_CURR=moist(IMS,KMS,JMS,P_QC) , F_QC=F_QC     &
1124      &               ,QR_CURR=moist(IMS,KMS,JMS,P_QR) , F_QR=F_QR     &
1125      &               ,QI_CURR=moist(IMS,KMS,JMS,P_QI) , F_QI=F_QI     &
1126      &               ,QS_CURR=moist(IMS,KMS,JMS,P_QS) , F_QS=F_QS     &
1127      &               ,QG_CURR=moist(IMS,KMS,JMS,P_QG) , F_QG=F_QG   )
1128 !
1129 !***  NOTE THAT THE EXCHANGE COEFFICIENTS FOR HEAT EXCH_H COMING OUT OF
1130 !***  PBL_DRIVER ARE DEFINED AT THE TOPS OF THE LAYERS KTS TO KTE-1
1131 !***  IF MODULE_BL_MYJPBL WAS INVOKED.
1132 !
1133 !-----------------------------------------------------------------------
1134 ! UNCOMPUTED LOCATIONS MUST BE FILLED IN FOR THE POST-PROCESSOR
1135 !-----------------------------------------------------------------------
1136 !
1137 !***  EASTERN GLOBAL BOUNDARY
1138 !
1139       IF(MYIE==IDE)THEN
1140 !$omp parallel do                                                       &
1141 !$omp& private(i,j)
1142         DO J=JDS,JDE
1143         IF (J>=MYJS.AND.J<=MYJE)THEN
1144           TH10(MYIE,J)=TH10(MYIE-1,J)
1145           Q10(MYIE,J)=Q10(MYIE-1,J)
1146           U10(MYIE,J)=U10(MYIE-1,J)
1147           V10(MYIE,J)=V10(MYIE-1,J)
1148           TSHLTR(MYIE,J)=TSHLTR(MYIE-1,J)
1149           QSHLTR(MYIE,J)=QSHLTR(MYIE-1,J)
1150         ENDIF
1151         ENDDO
1152       ENDIF
1153 !
1154 !***  SOUTHERN GLOBAL BOUNDARY
1155 !
1156 
1157       IF(MYJS==1)THEN
1158         DO J=1,2
1159         DO I=IDS,IDE
1160           IF (I>=MYIS.AND.I<=MYIE) THEN
1161             TH10(I,J)=TH10(I,MYJS+2)
1162             Q10(I,J)=Q10(I,MYJS+2)
1163             U10(I,J)=U10(I,MYJS+2)
1164             V10(I,J)=V10(I,MYJS+2)
1165             TSHLTR(I,J)=TSHLTR(I,MYJS+2)
1166             QSHLTR(I,J)=QSHLTR(I,MYJS+2)
1167           ENDIF
1168         ENDDO
1169         ENDDO
1170       ENDIF
1171 !
1172 !***  NORTHERN GLOBAL BOUNDARY
1173 !
1174       IF(MYJE==JDE)THEN
1175 !$omp parallel do                                                       &
1176 !$omp& private(i,j)
1177         DO J=MYJE-1,MYJE
1178         DO I=IDS,IDE
1179           IF (I>=MYIS.AND.I<=MYIE) THEN
1180             TH10(I,J)=TH10(I,MYJE-2)
1181             Q10(I,J)=Q10(I,MYJE-2)
1182             U10(I,J)=U10(I,MYJE-2)
1183             V10(I,J)=V10(I,MYJE-2)
1184             TSHLTR(I,J)=TSHLTR(I,MYJE-2)
1185             QSHLTR(I,J)=QSHLTR(I,MYJE-2)
1186           ENDIF
1187         ENDDO
1188         ENDDO
1189       ENDIF
1190 !
1191       IF(CONFIG_FLAGS%SF_SFCLAY_PHYSICS==1)THEN ! non-NMM package
1192 !$omp parallel do                                                       &
1193 !$omp& private(i,j)
1194         DO J=MYJS1,MYJE1
1195         DO I=MYIS,MYIE1
1196 !         TSHLTR(I,J)=TSHLTR(I,J)*(1.E5/PSHLTR(I,J))**RCP
1197           IF(TSHLTR(I,J)<200..OR.TSHLTR(I,J)>350.)THEN
1198             WRITE(0,*)'Troublesome TSHLTR...I,J,TSHLTR,PSHLTR: ',       &
1199                I,J,TSHLTR(I,J),PSHLTR(I,J)
1200           ENDIF
1201 	ENDDO
1202 	ENDDO
1203       ENDIF
1204 !
1205 !-----------------------------------------------------------------------
1206 !***  COMPUTE MODEL LAYER CONTAINING THE TOP OF THE BOUNDARY LAYER
1207 !-----------------------------------------------------------------------
1208 !
1209       IF(CONFIG_FLAGS%BL_PBL_PHYSICS/=MYJPBLSCHEME)THEN
1210         LENGTH_ROW=MYIE1-MYIS1+1
1211         DO J=MYJS2,MYJE2
1212         DO I=MYIS1,MYIE1
1213           KPBL(I,J)=-1000
1214         ENDDO
1215         ENDDO
1216 !
1217 !$omp parallel do                                                       &
1218 !$omp& private(altitude,i,j,k,kount_all)
1219         DO J=MYJS2,MYJE2
1220           KOUNT_ALL=0
1221           find_kpbl : DO K=KTS,KTE
1222           DO I=MYIS1,MYIE1
1223             ALTITUDE=Z(I,K+1,J)-SFCZ(I,J)
1224             IF(PBLH(I,J)<=ALTITUDE.AND.KPBL(I,J)<0)THEN
1225               KPBL(I,J)=K
1226               KOUNT_ALL=KOUNT_ALL+1
1227             ENDIF
1228             IF(KOUNT_ALL==LENGTH_ROW)EXIT find_kpbl
1229           ENDDO
1230           ENDDO find_kpbl
1231         ENDDO
1232       ENDIF
1233 !
1234       IF(MODEL_CONFIG_REC%SF_SURFACE_PHYSICS(GRID%ID)==99)THEN
1235         SNO_FACTR=1.
1236       ELSE
1237         SNO_FACTR=1000.
1238       ENDIF
1239 !
1240 !$omp parallel do                                                       &
1241 !$omp& private(i,j)
1242       DO J=MYJS2,MYJE2
1243       DO I=MYIS1,MYIE1
1244         SNO(I,J)=SNOW(I,J)
1245         SI(I,J)=SNOWH(I,J)*SNO_FACTR
1246         LPBL(I,J)=KTE-KPBL(I,J)+1
1247       ENDDO
1248       ENDDO
1249 !
1250 !-----------------------------------------------------------------------
1251 !***  DIAGNOSTIC RADIATION ACCUMULATION
1252 !-----------------------------------------------------------------------
1253 !
1254 !$omp parallel do                                                       &
1255 !$omp& private(i,j,tsfc2)
1256       DO J=MYJS2,MYJE2
1257       DO I=MYIS,MYIE
1258         ASWIN (I,J)=ASWIN (I,J)+RSWIN(I,J)*HBM2(I,J)*FACTRS(I,J)
1259         ASWOUT(I,J)=ASWOUT(I,J)-RSWOUT(I,J)*HBM2(I,J)*FACTRS(I,J)
1260         ASWTOA(I,J)=ASWTOA(I,J)+RSWTOA(I,J)*HBM2(I,J)*FACTRS(I,J)
1261         ALWIN (I,J)=ALWIN (I,J)+RLW_DN_SFC(I,J)
1262         ALWOUT(I,J)=ALWOUT(I,J)-RADOT (I,J)*HBM2(I,J)
1263         ALWTOA(I,J)=ALWTOA(I,J)+RLWTOA(I,J)*HBM2(I,J)
1264 !
1265         TSFC2=TSFC(I,J)*TSFC(I,J)
1266         RADOT(I,J)=HBM2(I,J)*EPSR(I,J)*STBOLT*TSFC2*TSFC2
1267         THS(I,J)=TSFC(I,J)*EXNSFC(I,J)
1268         PREC(I,J)=0.
1269       ENDDO
1270       ENDDO
1271 !
1272 !-----------------------------------------------------------------------
1273 !***  UPDATE TEMPERATURE, SPECIFIC HUMIDITY, CLOUD, AND TKE.
1274 !-----------------------------------------------------------------------
1275 !
1276       E_BDY=(ITE>=IDE)
1277 !
1278 !$omp parallel do                                                       &
1279 !$omp& private(dqdt,dtdt,i,iend,j,k,qi,qold,qr,qw,ratiomx,i_m)
1280       DO J=MYJS2,MYJE2
1281         IEND=MYIE1
1282         IF(E_BDY.AND.MOD(J,2)==0)IEND=IEND-1
1283 !
1284         DO K=KTS,KTE
1285         DO I=MYIS1,IEND
1286           DTDT=RTHBLTEN(I,K,J)*PI_PHY(I,K,J)
1287           DQDT=RQVBLTEN(I,K,J)         !Mixing ratio tendency
1288           T(I,K,J)=T(I,K,J)+DTDT*DTPHS
1289           QOLD=Q(I,K,J)
1290           RATIOMX=QOLD/(1.-QOLD)+DQDT*DTPHS
1291           Q(I,K,J)=RATIOMX/(1.+RATIOMX)
1292 !         Q(I,K,J)=MAX(Q(I,K,J),EPSQ)
1293           QW=max(0.,MOIST(I,K,J,P_QC)+RQCBLTEN(I,K,J)*DTPHS )
1294           QI=max(0.,MOIST(I,K,J,P_QS)+RQIBLTEN(I,K,J)*DTPHS )
1295           QR=max(0.,MOIST(I,K,J,P_QR) )
1296 !          CWM(I,K,J)=QW+QI+QR
1297           CWM(I,K,J)=0. 
1298 !
1299           DO I_M=1,N_MOIST
1300              IF(I_M/=P_QV)THEN
1301                CWM(I,K,J)=CWM(I,K,J)+MOIST(I,K,J,I_M)
1302              ENDIF
1303              IF(I_M==P_QV)THEN
1304                MOIST(I,K,J,P_QV)=MAX(EPSQ,(MOIST(I,K,J,P_QV) + RQVBLTEN(I,K,J)*DTPHS) )
1305              ELSEIF (I_M==P_QC)THEN
1306                 CWM(I,K,J)=MAX(0., (CWM(I,K,J) + RQCBLTEN(I,K,J)*DTPHS) )
1307              ELSEIF(I_M==P_QI)THEN
1308                 CWM(I,K,J)=MAX(0., (CWM(I,K,J) + RQIBLTEN(I,K,J)*DTPHS) )
1309              ENDIF
1310           ENDDO
1311 !
1312           MOIST(I,K,J,P_QC)=QW
1313           MOIST(I,K,J,P_QS)=QI
1314           MOIST(I,K,J,P_QR)=QR
1315 !
1316           IF(CONFIG_FLAGS%MP_PHYSICS==ETAMPNEW)THEN
1317             IF(QI<=EPSQ)THEN  
1318               F_ICE(I,K,J)=0.
1319             ELSE
1320               F_ICE(I,K,J)=MAX(0.,MIN(1.,QI/CWM(I,K,J)))
1321             ENDIF
1322 !
1323             IF(QR<=EPSQ)THEN
1324               F_RAIN(I,K,J)=0.
1325             ELSE
1326               F_RAIN(I,K,J)=QR/(QW+QR)
1327             ENDIF
1328           ENDIF
1329 !
1330           Q2(I,K,J)=2.*TKE(I,K,J)
1331         ENDDO
1332         ENDDO
1333 !
1334       ENDDO
1335 !
1336 !-----------------------------------------------------------------------
1337 !***
1338 !***  SAVE SURFACE-RELATED FIELDS.
1339 !***
1340 !-----------------------------------------------------------------------
1341 !$omp parallel do                                                       &
1342 !$omp& private(i,j,llij,xlvrw)
1343       DO J=MYJS2,MYJE2
1344       DO I=MYIS1,MYIE1
1345         LLIJ=LOWLYR(I,J)
1346 !
1347 !-----------------------------------------------------------------------
1348 !***  INSTANTANEOUS SENSIBLE AND LATENT HEAT FLUX
1349 !-----------------------------------------------------------------------
1350 !
1351         TWBS(I,J)=-TWBS(I,J)
1352         QWBS(I,J)=-QWBS(I,J)*XLV*CHKLOWQ(I,J)
1353 !
1354 !-----------------------------------------------------------------------
1355 !***  ACCUMULATED QUANTITIES.
1356 !***  IN OPNL LSM, SFCEVP APPEARS TO BE IN UNITS OF
1357 !***  METERS OF LIQUID WATER.  IT IS COMING FROM
1358 !***  WRF MODULE AS KG/M**2.
1359 !-----------------------------------------------------------------------
1360 !
1361         SFCSHX(I,J)=SFCSHX(I,J)+TWBS(I,J)
1362         SFCLHX(I,J)=SFCLHX(I,J)+QWBS(I,J)
1363         XLVRW=DTPHS/(XLV*RHOWATER)
1364         SFCEVP(I,J)=SFCEVP(I,J)-QWBS(I,J)*XLVRW
1365         POTEVP(I,J)=POTEVP(I,J)-QWBS(I,J)*SM(I,J)*XLVRW
1366         POTFLX(I,J)=POTEVP(I,J)*FACTOR
1367         SUBSHX(I,J)=SUBSHX(I,J)+GRNFLX(I,J)
1368       ENDDO
1369       ENDDO
1370 !
1371 !-----------------------------------------------------------------------
1372 !***  COUNTERS
1373 !-----------------------------------------------------------------------
1374 !
1375       APHTIM=APHTIM+1.
1376       ARDSW =ARDSW +1.
1377       ARDLW =ARDLW +1.
1378       ASRFC =ASRFC +1.
1379 !-----------------------------------------------------------------------
1380 !
1381       END SUBROUTINE TURBL
1382 !
1383 !-----------------------------------------------------------------------
1384 !***********************************************************************
1385       SUBROUTINE UV_H_TO_V(NTSD,DT,NPHS,UZ0H,VZ0H,UZ0,VZ0               &
1386      &                    ,DUDT,DVDT,U,V,HBM2,VTM,IVE,IVW               & 
1387      &                    ,IDS,IDE,JDS,JDE,KDS,KDE                      &
1388      &                    ,IMS,IME,JMS,JME,KMS,KME                      &
1389      &                    ,ITS,ITE,JTS,JTE,KTS,KTE)
1390 !***********************************************************************
1391 !$$$  SUBPROGRAM DOCUMENTATION BLOCK
1392 !                .      .    .     
1393 ! SUBPROGRAM:    UV_H_TO_V   INTERPOLATE WINDS FROM H TO V POINTS
1394 !   PRGRMMR: BLACK           ORG: W/NP22     DATE: 05-02-22       
1395 !     
1396 ! ABSTRACT:
1397 !     INTERPOLATE WINDS BACK TO V POINTS AFTER TURBULENCE
1398 !     
1399 ! PROGRAM HISTORY LOG :
1400 !   05-02-22  BLACK      - ORIGINATOR
1401 !     
1402 ! USAGE: CALL TURBL FROM SOLVE_NMM
1403 !
1404 ! ATTRIBUTES:
1405 !   LANGUAGE: FORTRAN 90
1406 !   MACHINE : IBM
1407 !$$$  
1408 !-----------------------------------------------------------------------
1409 !
1410       IMPLICIT NONE
1411 !
1412 !-----------------------------------------------------------------------
1413 !
1414       INTEGER,INTENT(IN) :: IDS,IDE,JDS,JDE,KDS,KDE                     &
1415      &                     ,IMS,IME,JMS,JME,KMS,KME                     &
1416      &                     ,ITS,ITE,JTS,JTE,KTS,KTE                     &
1417      &                     ,NPHS,NTSD
1418 !
1419       INTEGER, DIMENSION(JMS:JME),INTENT(IN) :: IVE,IVW
1420 !
1421       REAL,INTENT(IN) :: DT
1422 !
1423       REAL,DIMENSION(IMS:IME,JMS:JME),INTENT(IN) :: HBM2,UZ0H,VZ0H
1424 !
1425       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(IN) :: DUDT,DVDT   &
1426      &                                                     ,VTM
1427 !
1428       REAL,DIMENSION(IMS:IME,JMS:JME),INTENT(OUT) :: UZ0,VZ0
1429 !
1430       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(INOUT) :: U,V
1431 !
1432 !-----------------------------------------------------------------------
1433 !***
1434 !***  LOCAL VARIABLES
1435 !***
1436 !-----------------------------------------------------------------------
1437 !
1438       INTEGER :: I,IEND,J,K
1439 !
1440       REAL :: DTPHS
1441 !
1442       LOGICAL :: E_BDY
1443 !
1444 !-----------------------------------------------------------------------
1445 !-----------------------------------------------------------------------
1446 !
1447       DTPHS=NPHS*DT
1448       E_BDY=(ITE>=IDE)
1449 !
1450 !-----------------------------------------------------------------------
1451 !***  RECONSTRUCT UZ0 AND VZ0 ON VELOCITY POINTS.
1452 !-----------------------------------------------------------------------
1453 !
1454 !$omp parallel do                                                       &
1455 !$omp& private(i,j)
1456       DO J=MYJS2,MYJE2
1457       DO I=MYIS,MYIE
1458         UZ0(I,J)=(UZ0H(I+IVE(J),J)*HBM2(I+IVE(J),J)                     &
1459      &           +UZ0H(I+IVW(J),J)*HBM2(I+IVW(J),J)                     &
1460      &           +UZ0H(I,J+1)*HBM2(I,J+1)+UZ0H(I,J-1)*HBM2(I,J-1))*0.25
1461         VZ0(I,J)=(VZ0H(I+IVE(J),J)*HBM2(I+IVE(J),J)                     &
1462      &           +VZ0H(I+IVW(J),J)*HBM2(I+IVW(J),J)                     &
1463      &           +VZ0H(I,J+1)*HBM2(I,J+1)+VZ0H(I,J-1)*HBM2(I,J-1))*0.25
1464       ENDDO
1465       ENDDO
1466 !
1467 !-----------------------------------------------------------------------
1468 !***  INTERPOLATE WIND TENDENCIES TO VELOCITY POINTS AND UPDATE WINDS.
1469 !-----------------------------------------------------------------------
1470 !
1471 !$omp parallel do                                                       &
1472 !$omp& private(i,iend,j,k)
1473       DO J=MYJS2,MYJE2
1474         IEND=MYIE1
1475         IF(E_BDY.AND.MOD(J,2)==1)IEND=IEND-1
1476 !
1477         DO K=KTS,KTE
1478         DO I=MYIS1,IEND
1479           U(I,K,J)=(DUDT(I+IVE(J),K,J)+DUDT(I+IVW(J),K,J)               &
1480      &             +DUDT(I,K,J+1)+DUDT(I,K,J-1))*0.25*DTPHS             &
1481      &             *VTM(I,K,J)+U(I,K,J)
1482           V(I,K,J)=(DVDT(I+IVE(J),K,J)+DVDT(I+IVW(J),K,J)               &
1483      &             +DVDT(I,K,J+1)+DVDT(I,K,J-1))*0.25*DTPHS             &
1484      &             *VTM(I,K,J)+V(I,K,J)
1485         ENDDO
1486         ENDDO
1487       ENDDO
1488 !-----------------------------------------------------------------------
1489 !
1490       END SUBROUTINE UV_H_TO_V
1491 !
1492 !-----------------------------------------------------------------------
1493 !***********************************************************************
1494       SUBROUTINE CUCNVC(NTSD,DT,NCNVC,NRADS,NRADL                       &
1495      &                 ,GPS,RESTRT,HYDRO                                &
1496      &                 ,CLDEFI,LMH,N_MOIST,ENSDIM                       &
1497      &                 ,MOIST                                           &
1498      &                 ,DETA1,DETA2,AETA1,AETA2,ETA1,ETA2               &
1499      &                 ,F_ICE,F_RAIN                                    &
1500 !***  Changes for other cu-schemes, most for gd scheme
1501      &                 ,APR_GR,APR_W,APR_MC,TTEN,QTEN                   &
1502      &                 ,APR_ST,APR_AS,APR_CAPMA                         &
1503      &                 ,APR_CAPME          ,APR_CAPMI                   &
1504      &                 ,MASS_FLUX         ,XF_ENS                       &
1505      &                 ,PR_ENS,GSW                                      &
1506 #ifdef WRF_CHEM
1507      &                 ,GD_CLOUD,GD_CLOUD2                              &
1508 #endif
1509 !
1510      &                 ,PDTOP,PT,PD,RES,PINT,T,Q,CWM,TCUCN              &
1511      &                 ,OMGALF,U,V,VTM,WINT,Z,FIS,W0AVG                 &
1512      &                 ,PREC,ACPREC,CUPREC,CUPPT,CPRATE                 &
1513      &                 ,SM,HBM2,LPBL,CNVBOT,CNVTOP                      &
1514      &                 ,HTOP,HBOT,HTOPD,HBOTD,HTOPS,HBOTS               &
1515      &                 ,RTHBLTEN,RQVBLTEN,RTHRATEN                      & 
1516      &                 ,AVCNVC,ACUTIM,ZERO_3D,IHE,IHW                   &
1517      &                 ,GRID,CONFIG_FLAGS                               &
1518      &                 ,IDS,IDE,JDS,JDE,KDS,KDE                         &
1519      &                 ,IMS,IME,JMS,JME,KMS,KME                         &
1520      &                 ,ITS,ITE,JTS,JTE,KTS,KTE)
1521 !***********************************************************************
1522 !$$$  SUBPROGRAM DOCUMENTATION BLOCK
1523 !                .      .    .     
1524 ! SUBPROGRAM:    CUCNVC      CONVECTIVE PRECIPITATION OUTER DRIVER
1525 !   PRGRMMR: BLACK           ORG: W/NP22     DATE: 02-03-21       
1526 !     
1527 ! ABSTRACT:
1528 !     CUCVNC DRIVES THE WRF CONVECTION SCHEMES
1529 !     
1530 ! PROGRAM HISTORY LOG:
1531 !   02-03-21  BLACK      - ORIGINATOR
1532 !   04-11-18  BLACK      - THREADED
1533 !     
1534 ! USAGE: CALL CUCNVC FROM SOLVE_NMM
1535 !
1536 ! ATTRIBUTES:
1537 !   LANGUAGE: FORTRAN 90
1538 !   MACHINE : IBM 
1539 !$$$  
1540 !-----------------------------------------------------------------------
1541 !
1542       IMPLICIT NONE
1543 !
1544 !-----------------------------------------------------------------------
1545 !
1546       INTEGER,INTENT(IN) :: ENSDIM                                      &
1547      &                     ,IDS,IDE,JDS,JDE,KDS,KDE                     &
1548      &                     ,IMS,IME,JMS,JME,KMS,KME                     &
1549      &                     ,ITS,ITE,JTS,JTE,KTS,KTE                     &
1550      &                     ,N_MOIST,NCNVC,NTSD,NRADS,NRADL
1551 !
1552       INTEGER, DIMENSION(JMS:JME),INTENT(IN) :: IHE,IHW
1553 !
1554       INTEGER,DIMENSION(IMS:IME,JMS:JME),INTENT(IN) :: LMH,LPBL
1555 !
1556       REAL,INTENT(IN) :: DT,GPS,PDTOP,PT
1557 !
1558       REAL,INTENT(INOUT) :: ACUTIM,AVCNVC
1559 !
1560       REAL,DIMENSION(KMS:KME-1),INTENT(IN) :: AETA1,AETA2,DETA1,DETA2
1561       REAL,DIMENSION(KMS:KME  ),INTENT(IN) :: ETA1,ETA2
1562 !
1563       REAL,DIMENSION(IMS:IME,JMS:JME),INTENT(IN) :: FIS,HBM2,PD,RES,SM
1564 !
1565       REAL,DIMENSION(IMS:IME,JMS:JME),INTENT(INOUT) :: ACPREC,CLDEFI    &
1566      &                                                ,CNVBOT,CNVTOP    &
1567      &                                                ,CUPPT,CUPREC     &
1568      &                                                ,HBOT,HTOP        &
1569      &                                                ,HBOTD,HTOPD      &
1570      &                                                ,HBOTS,HTOPS      &
1571      &                 ,APR_GR,APR_W,APR_MC                             &
1572      &                 ,APR_ST,APR_AS,APR_CAPMA                         &
1573      &                 ,APR_CAPME          ,APR_CAPMI                   &
1574      &                 ,MASS_FLUX                                       &
1575      &                 ,GSW  ,PREC,CPRATE
1576 !
1577       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(IN) :: F_ICE       &
1578      &                                                     ,F_RAIN
1579       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(INOUT) :: TTEN     &
1580      &                                                     ,QTEN        &
1581      &                                       ,RTHBLTEN,RQVBLTEN,RTHRATEN
1582 
1583 !
1584       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(INOUT) :: Q,T      &
1585      &                                                        ,CWM      &
1586      &                                                        ,TCUCN    &
1587      &                                                        ,W0AVG    &
1588      &                                                        ,WINT
1589 !
1590       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(IN) :: OMGALF      &
1591      &                                                     ,PINT,U,V    &
1592      &                                                     ,VTM,Z
1593 !
1594       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(INOUT) :: ZERO_3D
1595       REAL,DIMENSION(IMS:IME,jMS:jME,1:ENSDIM),INTENT(INOUT) ::          &
1596      &                         XF_ENS                                    &
1597      &                        ,PR_ENS
1598       
1599 !    
1600       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME,N_MOIST)                   &
1601      &                                           ,INTENT(INOUT) :: moist
1602 #ifdef WRF_CHEM
1603       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(INOUT) :: GD_CLOUD,GD_CLOUD2
1604 #endif
1605 !
1606       LOGICAL,INTENT(IN) :: HYDRO,RESTRT
1607 !
1608       TYPE(DOMAIN),TARGET :: GRID
1609 !
1610       TYPE(GRID_CONFIG_REC_TYPE),INTENT(IN) :: CONFIG_FLAGS
1611 !
1612 !-----------------------------------------------------------------------
1613 !***
1614 !***  LOCAL VARIABLES
1615 !***
1616 !-----------------------------------------------------------------------
1617       INTEGER :: I,ICLDCK,IENDX,J,K,MNTO,NCUBOT,NCUTOP,NSTEP_CNV        &
1618      &          ,N_TIMSTPS_OUTPUT
1619 !
1620       INTEGER,DIMENSION(IMS:IME,JMS:JME) :: KPBL,LBOT,LOWLYR,LTOP
1621 !
1622       REAL :: CAPA,CF_HI,DPL,DQDT,DTCNVC,DTDT,FICE,FRAIN,G_INV          &
1623      &       ,PCPCOL,PDSL,PLYR,QI,QL_K,QR,QW,RDTCNVC,RWMSK,WMSK,WC
1624 !
1625       REAL,DIMENSION(KMS:KME-1) :: QL,TL
1626 !
1627       REAL,DIMENSION(IMS:IME,JMS:JME) :: CUBOT,CUTOP,NCA,RAINC,RAINCV   &
1628      &                                  ,SFCZ,XLAND
1629 !
1630       REAL,DIMENSION(IMS:IME,KMS:KME) :: WMID
1631 !
1632       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME) :: DZ,P8W,P_PHY,PI_PHY    &
1633      &                                          ,RQCCUTEN,RQRCUTEN      &
1634      &                                          ,RQICUTEN,RQSCUTEN      &
1635      &                                          ,RQVCUTEN,RR,RTHCUTEN   &
1636      &                                          ,T_PHY,TH_PHY           &
1637      &                                          ,U_PHY,V_PHY
1638 !
1639       REAL,DIMENSION(IMS:IME,JMS:JME)        :: ZERO_2D
1640       REAL,DIMENSION(IMS:IME,JMS:JME,ENSDIM) :: ZERO_GD
1641 !
1642       LOGICAL :: RESTART,WARM_RAIN
1643       LOGICAL,DIMENSION(IMS:IME,JMS:JME) :: CU_ACT_FLAG
1644 !
1645 !-----------------------------------------------------------------------
1646 !***  FOR TEMPERATURE CHANGE CHECK ONLY.
1647 !-----------------------------------------------------------------------
1648       INTEGER :: DTEMP_CHECK=1.0
1649       REAL :: TCHANGE
1650 !-----------------------------------------------------------------------
1651 !***********************************************************************
1652 !-----------------------------------------------------------------------
1653 !***  RESET THE HBOT/HTOP CONVECTIVE CLOUD BOTTOM (BASE) AND TOP ARRAYS
1654 !***  USED IN RADIATION.  THEY STORE THE MAXIMUM VERTICAL LIMITS OF 
1655 !***  CONVECTIVE CLOUD BETWEEN RADIATION CALLS.  CUPPT IS THE ACCUMULATED
1656 !***  CONVECTIVE PRECIPITATION BETWEEN RADIATION CALLS.
1657 !-----------------------------------------------------------------------
1658 !
1659       IF(MOD(NTSD,NRADS)==0.OR.MOD(NTSD,NRADL)==0)THEN
1660          DO J=JMS,JME
1661          DO I=IMS,IME
1662            HTOP(I,J)=0.
1663            HBOT(I,J)=REAL(KTE+1)
1664            CUPPT(I,J)=0.
1665          ENDDO
1666          ENDDO
1667       ENDIF
1668 !-----------------------------------------------------------------------
1669       IF(MOD(NTSD,NCNVC)/=0.AND.                                      &
1670      &   CONFIG_FLAGS%CU_PHYSICS==BMJSCHEME)RETURN
1671       IF(MOD(NTSD,NCNVC)/=0.AND.                                      &
1672      &   CONFIG_FLAGS%CU_PHYSICS==SASSCHEME)RETURN
1673 !-----------------------------------------------------------------------
1674       NSTEP_CNV=NCNVC
1675 !
1676       RESTART=RESTRT
1677 !-----------------------------------------------------------------------
1678       IF(CONFIG_FLAGS%CU_PHYSICS==KFETASCHEME)THEN
1679 !
1680         IF(.NOT.RESTART.AND.NTSD==0)THEN
1681 !$omp parallel do                                                       &
1682 !$omp& private(i,j,k)
1683           DO J=JTS,JTE
1684           DO K=KTS,KTE
1685           DO I=ITS,ITE
1686             W0AVG(I,K,J)=0.
1687           ENDDO
1688           ENDDO
1689           ENDDO
1690         ENDIF
1691 !
1692       ENDIF
1693 !
1694 !-----------------------------------------------------------------------
1695 !***  GENERAL PREPARATION 
1696 !-----------------------------------------------------------------------
1697 !
1698       AVCNVC=AVCNVC+1.
1699       ACUTIM=ACUTIM+1.
1700 !
1701       DTCNVC=NCNVC*DT
1702       RDTCNVC=1./DTCNVC
1703       CAPA=R_D/CP
1704       G_INV=1./G
1705 !
1706 !$omp parallel do                                                       &
1707 !$omp& private(dpl,fice,frain,i,j,k,pdsl,plyr,ql,tl)
1708       DO J=MYJS2,MYJE2
1709       DO I=MYIS1,MYIE1
1710 !
1711         PDSL=PD(I,J)*RES(I,J)
1712         RAINCV(I,J)=0.
1713         RAINC(I,J)=0.
1714         P8W(I,KTS,J)=PD(I,J)+PDTOP+PT
1715         LOWLYR(I,J)=KTE+1-LMH(I,J)
1716         XLAND(I,J)=SM(I,J)+1.
1717         NCA(I,J)=0.
1718         SFCZ(I,J)=FIS(I,J)*G_INV
1719 !tgs
1720           CUTOP(I,J)=HTOP(I,J)
1721           CUBOT(I,J)=HBOT(I,J)
1722 !
1723 !***  LPBL IS THE MODEL LAYER CONTAINING THE PBL TOP
1724 !***  COUNTING DOWNWARD FROM THE TOP OF THE DOMAIN
1725 !***  SO KPBL IS THE SAME LAYER COUNTING UPWARD FROM 
1726 !***  THE GROUND.
1727 !
1728         KPBL(I,J)=KTE-LPBL(I,J)+1
1729         ZERO_2D(I,J)=0
1730 !
1731         DO K=KTS,KTE
1732           DPL=DETA1(K)*PDTOP+DETA2(K)*PDSL
1733           QL(K)=AMAX1(Q(I,K,J),EPSQ)
1734           PLYR=AETA1(K)*PDTOP+AETA2(K)*PDSL+PT
1735           TL(K)=T(I,K,J)
1736 !
1737           RR(I,K,J)=PLYR/(R_D*TL(K)*(P608*QL(K)+1.))
1738           T_PHY(I,K,J)=TL(K)
1739 
1740           TH_PHY(I,K,J)=TL(K)*(1.E5/PLYR)**CAPA
1741 !!!       P8W(I,KFLIP,J)=PINT(I,K+1,J)
1742           P8W(I,K+1,J)=ETA1(K+1)*PDTOP+ETA2(K+1)*PDSL+PT
1743           P_PHY(I,K,J)=PLYR
1744           PI_PHY(I,K,J)=(PLYR*1.E-5)**CAPA
1745 !
1746           RTHCUTEN(I,K,J)=0.
1747           RQVCUTEN(I,K,J)=0.
1748           RQCCUTEN(I,K,J)=0.
1749           RQRCUTEN(I,K,J)=0.
1750           RQICUTEN(I,K,J)=0.
1751           RQSCUTEN(I,K,J)=0.
1752         ENDDO
1753 !
1754       ENDDO
1755       ENDDO
1756 !
1757 !-----------------------------------------------------------------------
1758 !
1759 
1760       IF(.NOT.HYDRO)THEN
1761 !$omp parallel do                                                       &
1762 !$omp& private(i,j,k)
1763         DO J=MYJS2,MYJE2
1764         DO K=KTS,KTE
1765         DO I=MYIS1,MYIE1
1766           DZ(I,K,J)=Z(I,K+1,J)-Z(I,K,J)
1767         ENDDO
1768         ENDDO
1769         ENDDO
1770 !
1771         IF(NTSD==0)THEN
1772 !$omp parallel do                                                       &
1773 !$omp& private(i,j,k)
1774           DO J=MYJS2,MYJE2
1775           DO K=KTS,KTE
1776           DO I=MYIS1,MYIE1
1777             WINT(I,K,J)=0.
1778           ENDDO
1779           ENDDO
1780           ENDDO
1781         ENDIF
1782       ELSE
1783         DO J=MYJS2,MYJE2
1784         DO I=MYIS1,MYIE1
1785           WINT(I,1,J)=0.
1786           WINT(I,KTE+1,J)=0.
1787         ENDDO
1788         ENDDO
1789 !
1790 !$omp parallel do                                                       &
1791 !$omp& private(i,j,k,plyr,wmid)
1792         DO J=MYJS2,MYJE2
1793           DO I=MYIS1,MYIE1
1794             WMID(I,KTS)=-OMGALF(I,KTS,J)*CP/(G*DT)
1795             PDSL=PD(I,J)*RES(I,J)
1796             PLYR=AETA1(KTS)*PDTOP+AETA2(KTS)*PDSL+PT
1797             DZ(I,KTS,J)=T(I,KTS,J)*(P608*Q(I,KTS,J)+1.)*R_D             &
1798      &                 *(P8W(I,KTS,J)-P8W(I,KTS+1,J))                   &
1799      &                 /(PLYR*G)
1800           ENDDO
1801 !
1802           DO K=KTS+1,KTE
1803           DO I=MYIS1,MYIE1
1804             QL_K=AMAX1(Q(I,K,J),EPSQ)
1805             WMID(I,K)=-OMGALF(I,K,J)*CP/(G*DT)
1806             WINT(I,K,J)=0.5*(WMID(I,K-1)+WMID(I,K))
1807             DZ(I,K,J)=T_PHY(I,K,J)*(P608*QL_K+1.)*R_D                   &
1808      &               *(P8W(I,K,J)-P8W(I,K+1,J))                         &
1809      &               /(P_PHY(I,K,J)*G)
1810           ENDDO
1811           ENDDO
1812         ENDDO
1813 !
1814       ENDIF
1815 !
1816 !-----------------------------------------------------------------------
1817 !***  COMPUTE VELOCITY COMPONENTS AT MASS POINTS
1818 !-----------------------------------------------------------------------
1819 !
1820       IF(CONFIG_FLAGS%CU_PHYSICS.NE.BMJSCHEME)THEN
1821 !
1822 !$omp parallel do                                                       &
1823 !$omp& private(i,j,k,rwmsk,wmsk)
1824         DO J=MYJS1_P1,MYJE1_P1
1825 !
1826           DO K=KTS,KTE
1827           DO I=MYIS_P1,MYIE_P1
1828             WMSK=VTM(I+IHE(J),K,J)+VTM(I+IHW(J),K,J)                    &
1829      &          +VTM(I,K,J+1)+VTM(I,K,J-1)
1830             IF(WMSK>0.)THEN
1831               RWMSK=1./WMSK
1832               U_PHY(I,K,J)=(U(I+IHE(J),K,J)*VTM(I+IHE(J),K,J)           &
1833      &                         +U(I+IHW(J),K,J)*VTM(I+IHW(J),K,J)       &
1834      &                         +U(I,K,J+1)*VTM(I,K,J+1)                 &
1835      &                         +U(I,K,J-1)*VTM(I,K,J-1))*RWMSK
1836               V_PHY(I,K,J)=(V(I+IHE(J),K,J)*VTM(I+IHE(J),K,J)           &
1837      &                         +V(I+IHW(J),K,J)*VTM(I+IHW(J),K,J)       &
1838      &                         +V(I,K,J+1)*VTM(I,K,J+1)                 &
1839      &                         +V(I,K,J-1)*VTM(I,K,J-1))*RWMSK
1840             ELSE
1841               U_PHY(I,K,J)=0.
1842               V_PHY(I,K,J)=0.
1843             ENDIF
1844           ENDDO
1845           ENDDO
1846 !
1847         ENDDO
1848 !
1849       ENDIF
1850 !-----------------------------------------------------------------------
1851 !
1852 !***  SINGLE-COLUMN CONVECTION
1853 !
1854 !-----------------------------------------------------------------------
1855 !
1856       CALL SET_TILES(GRID,IDS+1,IDE-1,JDS+2,JDE-2,ITS,ITE,JTS,JTE)
1857 !
1858       CALL CUMULUS_DRIVER(                                              &
1859      &                  IDS=IDS,IDE=IDE,JDS=JDS,JDE=JDE,KDS=KDS,KDE=KDE &
1860      &                 ,IMS=IMS,IME=IME,JMS=JMS,JME=JME,KMS=KMS,KME=KME &
1861      &                 ,I_START=GRID%I_START,I_END=GRID%I_END           &
1862      &                 ,J_START=GRID%J_START,J_END=GRID%J_END           &
1863      &                 ,KTS=KTS,KTE=KTE,NUM_TILES=GRID%NUM_TILES        &
1864                   ! Prognostic
1865      &                 ,U=U_PHY,V=V_PHY,TH=TH_PHY,T=T_PHY,W=WINT        &
1866      &                 ,P=P_PHY,PI=PI_PHY,RHO=RR,W0AVG=W0AVG            &
1867                   ! Others
1868      &                 ,ITIMESTEP=NTSD,DT=DT,DX=GPS                     &
1869      &                 ,RAINC=RAINC,RAINCV=RAINCV,NCA=NCA               &
1870      &                 ,DZ8W=DZ,P8W=P8W,FORCET=TTEN,FORCEQ=QTEN         &
1871      &                 ,CLDEFI=cldefi,LOWLYR=lowlyr,XLAND=xland         &
1872      &                 ,CU_ACT_FLAG=cu_act_flag,WARM_RAIN=warm_rain     &
1873      &                 ,STEPCU=NSTEP_CNV,GSW=gsw                        &
1874      &                 ,HTOP=CUTOP,HBOT=CUBOT,KPBL=KPBL,HT=SFCZ         &   
1875      &                 ,APR_GR=apr_gr,APR_W=apr_w,APR_MC=apr_mc         &
1876      &                 ,APR_ST=apr_st,APR_AS=apr_as,APR_CAPMA=apr_capma &
1877      &                 ,APR_CAPME=apr_capme,APR_CAPMI=apr_capmi         &
1878      &                 ,MASS_FLUX=mass_flux,XF_ENS=xf_ens               &
1879      &                 ,PR_ENS=pr_ens                                   &
1880 #ifdef WRF_CHEM
1881      &                 ,gd_cloud=gd_cloud,gd_cloud2=gd_cloud2           &
1882 #endif
1883 
1884      &                 ,ENSDIM=ENSDIM,MAXIENS=1,MAXENS=3                &
1885      &                 ,MAXENS2=3,MAXENS3=16                            &
1886      &                 ,RTHCUTEN=RTHCUTEN ,RQVCUTEN=RQVCUTEN            &
1887      &                 ,RQCCUTEN=RQCCUTEN ,RQRCUTEN=RQRCUTEN            &
1888      &                 ,RQICUTEN=RQICUTEN ,RQSCUTEN=RQSCUTEN            &
1889      &                 ,RTHBLTEN=RTHBLTEN,RQVBLTEN=RQVBLTEN             & 
1890      &                 ,RTHRATEN=RTHRATEN                               & 
1891                   ! Selection argument
1892      &                 ,CU_PHYSICS=CONFIG_FLAGS%CU_PHYSICS              &
1893                   ! Moisture tracer arguments
1894      &                 ,QV_CURR=MOIST(IMS,KMS,JMS,P_QV),F_QV=F_QV       &
1895      &                 ,QC_CURR=MOIST(IMS,KMS,JMS,P_QC),F_QC=F_QC       &
1896      &                 ,QR_CURR=MOIST(IMS,KMS,JMS,P_QR),F_QR=F_QR       &
1897      &                 ,QI_CURR=MOIST(IMS,KMS,JMS,P_QI),F_QI=F_QI       &
1898      &                 ,QS_CURR=MOIST(IMS,KMS,JMS,P_QS),F_QS=F_QS       &
1899      &                 ,QG_CURR=MOIST(IMS,KMS,JMS,P_QG),F_QG=F_QG      )
1900 !
1901 !-----------------------------------------------------------------------
1902 !
1903 !***  CNVTOP/CNVBOT HOLD THE MAXIMUM VERTICAL LIMITS OF CONVECTIVE CLOUD 
1904 !***  BETWEEN HISTORY OUTPUT TIMES.  HBOTS/HTOPS STORE SIMILIAR INFORMATION
1905 !***  FOR SHALLOW (NONPRECIPITATING) CONVECTION, AND HBOTD/HTOPD ARE FOR
1906 !***  DEEP (PRECIPITATING) CONVECTION.  
1907 !
1908       CF_HI=CONFIG_FLAGS%HISTORY_INTERVAL
1909       N_TIMSTPS_OUTPUT=NINT(60.*CF_HI/DT)
1910       MNTO=MOD(NTSD,N_TIMSTPS_OUTPUT)
1911 !
1912       IF(MNTO>0.AND.MNTO<=NCNVC)THEN
1913         DO J=MYJS2,MYJE2
1914         IENDX=MYIE1
1915         IF(MOD(J,2)==0.AND.ITE==IDE-1)IENDX=IENDX-1
1916         DO I=MYIS1,IENDX
1917           CNVBOT(I,J)=REAL(KTE+1.)
1918           CNVTOP(I,J)=0.
1919           HBOTD(I,J)=REAL(KTE+1.)
1920           HTOPD(I,J)=0.
1921           HBOTS(I,J)=REAL(KTE+1.)
1922           HTOPS(I,J)=0.
1923         ENDDO
1924         ENDDO
1925       ENDIF
1926 !
1927 !-----------------------------------------------------------------------
1928 !
1929 !$omp parallel do                                                       &
1930 !$omp& private(dqdt,dtdt,i,iendx,j,k,ncubot,ncutop,pcpcol               &
1931 !$omp&        ,tchange                                                  &
1932 !$omp&        )
1933       DO J=MYJS2,MYJE2
1934       IENDX=MYIE1
1935       IF(MOD(J,2)==0.AND.ITE==IDE-1)IENDX=IENDX-1
1936       DO I=MYIS1,IENDX
1937 !
1938 !***  UPDATE TEMPERATURE, SPECIFIC HUMIDITY, AND HEATING.
1939 !***  THE FLIP IS BECAUSE RTHCUTEN AND RQVCUTEN REACH THIS POINT
1940 !***  WITH LAYER 1 AT THE BOTTOM.
1941 !
1942         DO K=KTS,KTE
1943 !
1944 !***  RQVCUTEN IN BMJDRV IS THE MIXING RATIO TENDENCY,
1945 !***  SO RETRIEVE DQDT BY CONVERTING TO SPECIFIC HUMIDITY.
1946 !
1947           DQDT=RQVCUTEN(I,K,J)/(1.+MOIST(I,K,J,P_QV))**2
1948 !
1949 !***  RTHCUTEN IN BMJDRV IS DTDT OVER PI.
1950 !
1951           DTDT=RTHCUTEN(I,K,J)*PI_PHY(I,K,J)
1952           T(I,K,J)=T(I,K,J)+DTDT*DTCNVC
1953           Q(I,K,J)=Q(I,K,J)+DQDT*DTCNVC
1954           MOIST(I,K,J,P_QV)=Q(I,K,J)/(1.-Q(I,K,J))       !Convert to mixing ratio
1955 !tgs - added next two lines
1956           cps_select: SELECT CASE(config_flags%cu_physics)
1957 !
1958           CASE (KFSCHEME,KFETASCHEME,GDSCHEME,SASSCHEME)
1959              MOIST(I,K,J,P_QS)=MAX(0.,MOIST(I,K,J,P_QS)+RQICUTEN(I,K,J)*DTCNVC)
1960              MOIST(I,K,J,P_QC)=MAX(0.,MOIST(I,K,J,P_QC)+RQCCUTEN(I,K,J)*DTCNVC)
1961           END SELECT cps_select
1962 !
1963           TCUCN(I,K,J)=TCUCN(I,K,J)+DTDT
1964 !
1965           TCHANGE=DTDT*DTCNVC
1966 	  IF(ABS(TCHANGE)>DTEMP_CHECK)THEN
1967             WRITE(0,*)'BIG T CHANGE BY CONVECTION:  I,J,K,NTSD',TCHANGE,I,J,K,NTSD
1968 	  ENDIF
1969 !
1970         ENDDO
1971 !
1972 !***  UPDATE PRECIPITATION
1973 !
1974         PCPCOL=RAINCV(I,J)*1.E-3*NSTEP_CNV
1975         PREC(I,J)=PREC(I,J)+PCPCOL
1976         ACPREC(I,J)=ACPREC(I,J)+PCPCOL
1977         CUPREC(I,J)=CUPREC(I,J)+PCPCOL
1978         CUPPT(I,J)=CUPPT(I,J)+PCPCOL
1979         CPRATE(I,J)=PCPCOL
1980 !
1981 !***  SAVE CLOUD TOP AND BOTTOM FOR RADIATION (HTOP/HBOT) AND
1982 !***  FOR OUTPUT (CNVTOP/CNVBOT, HTOPS/HBOTS, HTOPD/HBOTD) ARRAYS.
1983 !***  MUST BE TREATED SEPARATELY FROM EACH OTHER.
1984 !
1985         NCUTOP=NINT(CUTOP(I,J))
1986         NCUBOT=NINT(CUBOT(I,J))
1987 !
1988         IF(NCUTOP>1.AND.NCUTOP<KDE)THEN
1989           HTOP(I,J)=MAX(CUTOP(I,J),HTOP(I,J))
1990           CNVTOP(I,J)=MAX(CUTOP(I,J),CNVTOP(I,J))
1991           IF(PCPCOL>0.)THEN
1992             HTOPD(I,J)=MAX(CUTOP(I,J),HTOPD(I,J))
1993           ELSE
1994             HTOPS(I,J)=MAX(CUTOP(I,J),HTOPS(I,J))
1995           ENDIF
1996         ENDIF
1997         IF(NCUBOT>0.AND.NCUBOT<KDE)THEN
1998           HBOT(I,J)=MIN(CUBOT(I,J),HBOT(I,J))
1999           CNVBOT(I,J)=MIN(CUBOT(I,J),CNVBOT(I,J))
2000           IF(PCPCOL>0.)THEN
2001             HBOTD(I,J)=MIN(CUBOT(I,J),HBOTD(I,J))
2002           ELSE
2003             HBOTS(I,J)=MIN(CUBOT(I,J),HBOTS(I,J))
2004           ENDIF
2005         ENDIF
2006 !
2007       ENDDO
2008       ENDDO
2009 !
2010 !$omp parallel do                                                       &
2011 !$omp& private(i,j,k)
2012       DO J=JMS,JME
2013       DO K=KMS,KME
2014       DO I=IMS,IME
2015         ZERO_3D(I,K,J)=0.
2016       ENDDO
2017       ENDDO
2018       ENDDO
2019 !-----------------------------------------------------------------------
2020 !
2021       END SUBROUTINE CUCNVC
2022 !
2023 !-----------------------------------------------------------------------
2024 !***********************************************************************
2025       SUBROUTINE GSMDRIVE(NTSD,DT,NPHS,N_MOIST                          &
2026      &                   ,DX,DY,LMH,SM,HBM2,FIS                         &
2027      &                   ,DETA1,DETA2,AETA1,AETA2,ETA1,ETA2             &
2028      &                   ,PDTOP,PT,PD,RES,PINT,T,Q,CWM,TRAIN            &
2029      &                   ,MOIST,SCALAR,N_SCALAR                         &
2030      &                   ,F_ICE,F_RAIN,F_RIMEF,SR                       &
2031      &                   ,PREC,ACPREC,AVRAIN,ZERO_3D                    &
2032      &                   ,MP_RESTART_STATE                              &
2033      &                   ,TBPVS_STATE                                   &
2034      &                   ,TBPVS0_STATE                                  &
2035      &                   ,GRID,CONFIG_FLAGS                             &
2036      &                   ,IDS,IDE,JDS,JDE,KDS,KDE                       &
2037      &                   ,IMS,IME,JMS,JME,KMS,KME                       &
2038      &                   ,ITS,ITE,JTS,JTE,KTS,KTE)
2039 !***********************************************************************
2040 !$$$  SUBPROGRAM DOCUMENTATION BLOCK
2041 !                .      .    .     
2042 ! SUBPROGRAM:    GSMDRIVE    MICROPHYSICS OUTER DRIVER
2043 !   PRGRMMR: BLACK           ORG: W/NP22     DATE: 02-03-26       
2044 !     
2045 ! ABSTRACT:
2046 !     GSMDRIVE DRIVES THE MICROPHYSICS SCHEMES
2047 !     
2048 ! PROGRAM HISTORY LOG:
2049 !   02-03-26  BLACK      - ORIGINATOR
2050 !   04-11-18  BLACK      - THREADED
2051 !     
2052 ! USAGE: CALL GSMDRIVE FROM SOLVE_NMM
2053 !
2054 ! ATTRIBUTES:
2055 !   LANGUAGE: FORTRAN 90
2056 !   MACHINE : IBM
2057 !$$$  
2058 !-----------------------------------------------------------------------
2059 !
2060       IMPLICIT NONE
2061 !
2062 !-----------------------------------------------------------------------
2063 !
2064       INTEGER,INTENT(IN) :: IDS,IDE,JDS,JDE,KDS,KDE                     &
2065      &                     ,IMS,IME,JMS,JME,KMS,KME                     &
2066      &                     ,ITS,ITE,JTS,JTE,KTS,KTE                     &
2067      &                     ,N_MOIST,N_SCALAR,NPHS,NTSD
2068 !
2069       INTEGER,DIMENSION(IMS:IME,JMS:JME),INTENT(IN) :: LMH
2070 !
2071       REAL,INTENT(IN) :: DT,DX,DY,PDTOP,PT
2072 !
2073       REAL,INTENT(INOUT) :: AVRAIN
2074 !
2075       REAL,DIMENSION(KMS:KME-1),INTENT(IN) :: AETA1,AETA2,DETA1,DETA2
2076       REAL,DIMENSION(KMS:KME),INTENT(IN) :: ETA1,ETA2
2077 !
2078       REAL,DIMENSION(IMS:IME,JMS:JME),INTENT(IN) :: FIS,HBM2,PD,RES,SM
2079 !
2080       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(IN) :: PINT
2081       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(INOUT) :: ZERO_3D
2082 !
2083       REAL,DIMENSION(IMS:IME,JMS:JME),INTENT(INOUT) :: ACPREC,PREC
2084 !
2085       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(INOUT) :: CWM,Q,T  &
2086      &                                                        ,TRAIN
2087 !
2088       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME),INTENT(INOUT) :: F_ICE    &
2089      &                                                        ,F_RAIN   &
2090      &                                                        ,F_RIMEF
2091 
2092       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME,n_moist),INTENT(INOUT) :: MOIST
2093       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME,n_scalar),INTENT(INOUT) :: SCALAR
2094 !
2095 !***  State var for etampnew microphysics (JM, 2005 05 02)
2096 !
2097       REAL,DIMENSION(:),INTENT(INOUT) ::               MP_RESTART_STATE &
2098      &                                                     ,TBPVS_STATE &
2099      &                                                    ,TBPVS0_STATE
2100 
2101 !
2102       REAL,DIMENSION(IMS:IME,JMS:JME),INTENT(OUT) :: SR
2103 !
2104       TYPE(DOMAIN),TARGET :: GRID
2105 !
2106       TYPE(GRID_CONFIG_REC_TYPE),INTENT(IN) :: CONFIG_FLAGS
2107 !
2108 !-----------------------------------------------------------------------
2109 !***
2110 !***  LOCAL VARIABLES
2111 !***
2112 !-----------------------------------------------------------------------
2113       INTEGER :: I,I_M,IENDX,J,K,IJ
2114 !
2115       INTEGER,DIMENSION(IMS:IME,JMS:JME) :: LOWLYR
2116 !
2117       REAL :: CAPA,DPL,DTPHS,PCPCOL,PDSL,PLYR,RDTPHS,RG,TNEW
2118 !
2119       REAL,DIMENSION(KMS:KME-1) :: QL,TL
2120 !
2121       REAL,DIMENSION(IMS:IME,JMS:JME) :: CUBOT,CUTOP,RAINNC,RAINNCV,XLAND      &
2122      &                                  ,ZERO_2D
2123 !
2124       REAL,DIMENSION(IMS:IME,KMS:KME,JMS:JME) :: DZ,P8W,P_PHY,PI_PHY    &
2125      &                                          ,RR,T_PHY,TH_PHY
2126 !
2127       LOGICAL :: E_BDY,F_QT,QT_PRESENT,WARM_RAIN
2128 !
2129 !-----------------------------------------------------------------------
2130 !***********************************************************************
2131 !-----------------------------------------------------------------------
2132 !
2133       IF(CONFIG_FLAGS%MP_PHYSICS==ETAMPNEW)THEN
2134         QT_PRESENT=.TRUE.
2135       ELSE
2136         QT_PRESENT=.FALSE.
2137       ENDIF
2138 !
2139       DTPHS=NPHS*DT
2140       RDTPHS=1./DTPHS
2141       CAPA=R_D/CP
2142       RG=1./G
2143       AVRAIN=AVRAIN+1.
2144 !
2145 !-----------------------------------------------------------------------
2146 !
2147 !***  PREPARE NEEDED ARRAYS
2148 !
2149 !-----------------------------------------------------------------------
2150 !$omp parallel do                                                       &
2151 !$omp& private(dpl,i,j,k,pdsl,plyr,ql,tl)
2152       DO J=MYJS2,MYJE2
2153       DO I=MYIS1,MYIE1
2154 !
2155         PDSL=PD(I,J)*RES(I,J)
2156         P8W(I,KTE+1,J)=PT
2157         LOWLYR(I,J)=KTE+1-LMH(I,J)
2158         XLAND(I,J)=SM(I,J)+1.
2159         ZERO_2D(I,J)=0.
2160 ! FILL RAINNC WITH ZERO (NORMALLY CONTAINS THE NONCONVECTIVE 
2161 !         ACCUMULATED RAIN BUT NOT YET USED BY NMM)
2162 ! COULD BE OBTAINED FROM ACPREC AND CUPREC (ACPREC-CUPREC) 
2163         RAINNC(I,J)=0.
2164 !
2165 !***  FILL THE SINGLE-COLUMN INPUT
2166 !
2167         DO K=KTS,KTE
2168           DPL=DETA1(K)*PDTOP+DETA2(K)*PDSL
2169           QL(K)=AMAX1(Q(I,K,J),EPSQ)
2170 !!!       PLYR=AETA1(K)*PDTOP+AETA2(K)*PDSL+PT
2171           PLYR=(PINT(I,K,J)+PINT(I,K+1,J))*0.5
2172           TL(K)=T(I,K,J)
2173 !
2174           RR(I,K,J)=PLYR/(R_D*TL(K)*(P608*QL(K)+1.))
2175           T_PHY(I,K,J)=TL(K)
2176           PI_PHY(I,K,J)=(PLYR*1.E-5)**CAPA
2177           TH_PHY(I,K,J)=TL(K)/PI_PHY(I,K,J)
2178 !!!       P8W(I,KFLIP,J)=PINT(I,K+1,J)
2179           P8W(I,K,J)=ETA1(K)*PDTOP+ETA2(K)*PDSL+PT
2180           P_PHY(I,K,J)=PLYR
2181           DZ(I,K,J)=DPL*RG/RR(I,K,J)
2182         ENDDO
2183 !
2184       ENDDO
2185       ENDDO
2186 !-----------------------------------------------------------------------
2187 !
2188 !***  CALL MICROPHYSICS
2189 !
2190 !-----------------------------------------------------------------------
2191 !
2192       CALL SET_TILES(GRID,IDS+1,IDE-1,JDS+2,JDE-2,ITS,ITE,JTS,JTE)
2193 !
2194       CALL MICROPHYSICS_DRIVER(                                         &
2195      &                  TH=TH_PHY,RHO=RR,PI_PHY=PI_PHY,P=P_PHY          &
2196      &                 ,RAINNC=RAINNC,RAINNCV=RAINNCV                   &
2197      &                 ,DZ8W=DZ,P8W=P8W,DT=DTPHS,DX=DX,DY=DY            &
2198      &                 ,MP_PHYSICS=CONFIG_FLAGS%MP_PHYSICS              &
2199      &                 ,SPECIFIED=CONFIG_FLAGS%SPECIFIED                &
2200      &                        .OR.CONFIG_FLAGS%NESTED                   &
2201      &                 ,SPEC_ZONE=0,WARM_RAIN=WARM_RAIN                 &
2202      &                 ,XLAND=XLAND,ITIMESTEP=NTSD-1                    &
2203      &                 ,F_ICE_PHY=F_ICE,F_RAIN_PHY=F_RAIN               &
2204      &                 ,F_RIMEF_PHY=F_RIMEF                             &
2205      &                 ,LOWLYR=LOWLYR,SR=SR                             &
2206      &                 ,QV_CURR=MOIST(IMS,KMS,JMS,P_QV),F_QV=F_QV       &
2207      &                 ,QC_CURR=MOIST(IMS,KMS,JMS,P_QC),F_QC=F_QC       &
2208      &                 ,QR_CURR=MOIST(IMS,KMS,JMS,P_QR),F_QR=F_QR       &
2209      &                 ,QI_CURR=MOIST(IMS,KMS,JMS,P_QI),F_QI=F_QI       &
2210      &                 ,QS_CURR=MOIST(IMS,KMS,JMS,P_QS),F_QS=F_QS       &
2211      &                 ,QG_CURR=MOIST(IMS,KMS,JMS,P_QG),F_QG=F_QG       &
2212      &                 ,QNI_CURR=SCALAR(IMS,KMS,JMS,P_QNI),F_QNI=F_QNI  &
2213      &                 ,QT_CURR=CWM,F_QT=qt_present                     &
2214      &                 ,MP_RESTART_STATE=MP_RESTART_STATE               &
2215      &                 ,TBPVS_STATE=TBPVS_STATE                         &
2216      &                 ,TBPVS0_STATE=TBPVS0_STATE                       &
2217      &                 ,IDS=IDS,IDE=IDE,JDS=JDS,JDE=JDE,KDS=KDS,KDE=KDE &
2218      &                 ,IMS=IMS,IME=IME,JMS=JMS,JME=JME,KMS=KMS,KME=KME &
2219      &                 ,I_START=GRID%I_START,I_END=GRID%I_END           &
2220      &                 ,J_START=GRID%J_START,J_END=GRID%J_END           &
2221      &                 ,KTS=KTS,KTE=KTE,NUM_TILES=GRID%NUM_TILES        &
2222                                                                         )
2223 
2224 !$omp parallel do                                                       &
2225 !$omp& private(ij)
2226       DO IJ=1,GRID%NUM_TILES
2227         CALL MICROPHYSICS_ZERO_OUT(                                     &
2228                      MOIST,N_MOIST,CONFIG_FLAGS                         &
2229                     ,IDS,IDE,JDS,JDE,KDS,KDE                            &
2230                     ,IMS,IME,JMS,JME,KMS,KME                            &
2231                     ,GRID%I_START(IJ),GRID%I_END(IJ)                    &
2232                     ,GRID%J_START(IJ),GRID%J_END(IJ)                    &
2233                     ,KTS,KTE                                       )
2234       ENDDO
2235 
2236 
2237 
2238 !
2239 !-----------------------------------------------------------------------
2240 !
2241       E_BDY=(ITE>=IDE)
2242 !
2243 !$omp parallel do                                                       &
2244 !$omp& private(i,iendx,j,k,pcpcol,tnew,i_m)
2245       DO J=MYJS2,MYJE2
2246       IENDX=MYIE1
2247       IF(E_BDY.AND.MOD(J,2)==0)IENDX=IENDX-1
2248       DO I=MYIS1,IENDX
2249 !
2250 !***  UPDATE TEMPERATURE, SPECIFIC HUMIDITY, CLOUD WATER, AND HEATING.
2251 !
2252         DO K=KTS,KTE
2253           TNEW=TH_PHY(I,K,J)*PI_PHY(I,K,J)
2254           TRAIN(I,K,J)=TRAIN(I,K,J)+(TNEW-T(I,K,J))*RDTPHS
2255           T(I,K,J)=TNEW
2256           Q(I,K,J)=MOIST(I,K,J,P_QV)/(1.+MOIST(I,K,J,P_QV)) !To s.h.
2257 !         CWM(I,K,J)=0.
2258 !         DO I_M=2,N_MOIST
2259 !           IF(I_M/=P_QV)THEN
2260 !             CWM(I,K,J)=CWM(I,K,J)+MOIST(I,K,J,I_M)
2261 !           ENDIF
2262 !         ENDDO
2263         ENDDO
2264 !
2265 !-----------------------------------------------------------------------
2266 !***  UPDATE PRECIPITATION
2267 !-----------------------------------------------------------------------
2268 !
2269         PCPCOL=RAINNCV(I,J)*1.E-3
2270         PREC(I,J)=PREC(I,J)+PCPCOL
2271         ACPREC(I,J)=ACPREC(I,J)+PCPCOL
2272 ! NOTE: RAINNC IS ACCUMULATED INSIDE MICROPHYSICS BUT NMM ZEROES IT OUT ABOVE
2273 !    SINCE IT IS ONLY A LOCAL ARRAY FOR NOW
2274 !
2275       ENDDO
2276       ENDDO
2277 !
2278 !-----------------------------------------------------------------------
2279 !$omp parallel do                                                       &
2280 !$omp& private(i,j,k)
2281       DO J=JMS,JME
2282       DO K=KMS,KME
2283       DO I=IMS,IME
2284         ZERO_3D(I,K,J)=0.
2285       ENDDO
2286       ENDDO
2287       ENDDO
2288 !-------------------------------------------------------------------
2289 !
2290       END SUBROUTINE GSMDRIVE
2291 !
2292 !-------------------------------------------------------------------
2293 !
2294       END MODULE MODULE_PHYSICS_CALLS
2295 !
2296 !-------------------------------------------------------------------