merge standard release WRF/WPS V3.0.1.1 into wrffire
[wrffire.git] / wrfv2_fire / dyn_nmm / solve_nmm.F
blobe56765f35ee95e9f8297ee10a5aed05e75d3db48
1 !-----------------------------------------------------------------------
3 !NCEP_MESO:MEDIATION_LAYER:SOLVER
5 !-----------------------------------------------------------------------
6 #include "nmm_loop_basemacros.h"
7 #include "nmm_loop_macros.h"
8 !-----------------------------------------------------------------------
10       SUBROUTINE SOLVE_NMM(GRID,CONFIG_FLAGS                            &
12 #include "dummy_args.inc"
14      &           )
15 !-----------------------------------------------------------------------
16       USE MODULE_DOMAIN,                ONLY : DOMAIN, GET_IJK_FROM_GRID
17       USE MODULE_CONFIGURE,             ONLY : GRID_CONFIG_REC_TYPE
18       USE MODULE_MODEL_CONSTANTS
19       USE MODULE_STATE_DESCRIPTION
20       USE MODULE_CTLBLK
21 !#ifdef DM_PARALLEL
22 !      USE MODULE_DM
23 !#endif
24       USE MODULE_IGWAVE_ADJUST,         ONLY: PDTE,PFDHT,DDAMP,VTOA
25       USE MODULE_ADVECTION,             ONLY: ADVE,VAD2,HAD2,VAD2_SCAL,HAD2_SCAL
26       USE MODULE_NONHY_DYNAM,           ONLY: EPS,VADZ,HADZ
27       USE MODULE_DIFFUSION_NMM,         ONLY: HDIFF
28       USE MODULE_BNDRY_COND,            ONLY: BOCOH,BOCOV
29       USE MODULE_PHYSICS_CALLS
30       USE MODULE_EXT_INTERNAL
31       USE MODULE_PRECIP_ADJUST
32       USE MODULE_NEST_UTIL     ! USEs module_MPP (contains MYPE,NPES,MPI_COMM_COMP)
33 #ifdef WRF_CHEM
34       USE MODULE_INPUT_CHEM_DATA,       ONLY: GET_LAST_GAS
35 #endif
36 !-----------------------------------------------------------------------
38       IMPLICIT NONE
40 !-----------------------------------------------------------------------
42 !***  INPUT DATA
44 !-----------------------------------------------------------------------
46       TYPE(DOMAIN),TARGET :: GRID
48 !***  DEFINITIONS OF DUMMY ARGUMENTS TO THIS ROUTINE (GENERATED FROM REGISTRY)
50 ! NOTE, REGISTRY NO LONGER GENERATES DUMMY ARGUMENTS OR DUMMY ARGUMENT
51 ! DECLARATIONS FOR RCONFIG ENTRIES. THEY ARE STILL PART OF STATE. ACCESS
52 ! TO THESE VARIABLES IS NOW THROUGH GRID STRUCTURE, AS MODIFIED BELOW.
53 ! AFFECTED VARIABLES: SIGMA, DT, NPHS, IDTAD, NRADS, NRADL, JULDAY,
54 ! JULYR, NUM_SOIL_LAYERS, NCNVC, ENSDIM, DY, AND SPEC_BDY_WIDTH.
55 ! JM, 20050819
57 !----------------------------
58 #include <dummy_decl.inc>
59 !----------------------------
61 !***  STRUCTURE THAT CONTAINS RUN-TIME CONFIGURATION (NAMELIST) DATA FOR DOMAIN
63       TYPE(GRID_CONFIG_REC_TYPE),INTENT(IN) :: CONFIG_FLAGS
64 #ifdef WRF_CHEM
65       INTEGER :: NUMGAS
66 #endif
68 !-----------------------------------------------------------------------
70 !***  LOCAL VARIABLES
72 !-----------------------------------------------------------------------
73       INTEGER :: IDS,IDE,JDS,JDE,KDS,KDE                                &
74      &          ,IMS,IME,JMS,JME,KMS,KME                                & 
75      &          ,IPS,IPE,JPS,JPE,KPS,KPE                                &
76      &          ,ITS,ITE,JTS,JTE,KTS,KTE
78       INTEGER :: I,ICLTEND,IDF,IRTN,J,JC,JDF,K,KDF,LB,N_MOIST &
79      &          ,NTSD_current,L
80       integer :: ierr
81       INTEGER,SAVE :: NTSD_restart
82 !     INTEGER :: MPI_COMM_COMP,MYPE,MYPROC,NPES
83       INTEGER :: MYPROC
84       INTEGER :: KVH,NTSD_rad,RC
85       INTEGER :: NUM_OZMIXM,NUM_AEROSOLC
87       REAL :: DT_INV,FICE,FRAIN,GPS,QI,QR,QW,WC
89       LOGICAL :: LAST_TIME,OPERATIONAL_PHYSICS
91       CHARACTER(80) :: MESSAGE
93 !***  For precip assimilation:
94       INTEGER :: ISTAT
95       REAL,ALLOCATABLE,SAVE,DIMENSION(:,:,:) :: PPTDAT
96 #ifdef WRF_CHEM
97       REAL,ALLOCATABLE,DIMENSION(:,:,:,:) :: CHEM_TRANS
98 #endif
100 !-----------------------------------------------------------------------
101 !***  For physics compatibility with other packages
102 !-----------------------------------------------------------------------
104       REAL,ALLOCATABLE,DIMENSION(:,:,:) :: TTEN,QTEN
105       REAL,ALLOCATABLE,DIMENSION(:,:,:) :: RTHRATEN,RTHBLTEN,RQVBLTEN
107 !-----------------------------------------------------------------------
109       LOGICAL wrf_dm_on_monitor
110       EXTERNAL wrf_dm_on_monitor
112 !-----------------------------------------------------------------------
113 !***  TIMING VARIABLES
114 !-----------------------------------------------------------------------
115       real,save :: solve_tim,exch_tim,pdte_tim,adve_tim,vtoa_tim        &
116      &,            vadz_tim,hadz_tim,eps_tim,vad2_tim,had2_tim          &
117      &,            radiation_tim,rdtemp_tim,turbl_tim,cltend_tim        &
118      &,            cucnvc_tim,gsmdrive_tim,hdiff_tim,bocoh_tim          &
119      &,            pfdht_tim,ddamp_tim,bocov_tim,uv_htov_tim,sum_tim    &
120      &,            adjppt_tim
121       real,save :: exch_tim_max
122       real :: btim,btimx
123       real :: et_max,this_tim
124       integer :: n_print_time
126 #ifdef RSL
127       integer rsl_internal_milliclock
128       external rsl_internal_milliclock
129 # define timef rsl_internal_milliclock
130 #else
131       real*8 :: timef
132 #endif
133 !-----------------------------------------------------------------------
135 #ifdef DEREF_KLUDGE
136 ! SEE http://www.mmm.ucar.edu/wrf/WG2/topics/deref_kludge.htm
137       INTEGER :: SM31,EM31,SM32,EM32,SM33,EM33 
138       INTEGER :: SM31X,EM31X,SM32X,EM32X,SM33X,EM33X
139       INTEGER :: SM31Y,EM31Y,SM32Y,EM32Y,SM33Y,EM33Y
140 #endif
142 !-----------------------------------------------------------------------
144 ! LIMIT THE NUMBER OF ARGUMENTS IF COMPILED WITH -DLIMIT_ARGS BY COPYING
145 ! SCALAR (NON-ARRAY) ARGUMENTS OUT OF THE GRID DATA STRUCTURE INTO LOCALLY
146 ! DEFINED COPIES (DEFINED IN EM_DUMMY_DECL.INC, ABOVE, AS THEY ARE IF THEY
147 ! ARE ARGUMENTS).  AN EQUIVALENT INCLUDE OF EM_SCALAR_DEREFS.INC APPEARS
148 ! AT THE END OF THE ROUTINE TO COPY BACK ANY CHNAGED NON-ARRAY VALUES.
149 ! THE DEFINITION OF COPY_IN OR COPY_OUT BEFORE THE INCLUDE DEFINES THE
150 ! DIRECTION OF THE COPY.  NMM_SCALAR_DEREFS IS GENERATED FROM REGISTRY.
152 !-----------------------------------------------------------------------
153 #define COPY_IN
154 #include <scalar_derefs.inc>
155 !-----------------------------------------------------------------------
157 ! TRICK PROBLEMATIC COMPILERS INTO NOT PERFORMING COPY-IN/COPY-OUT BY ADDING
158 ! INDICES TO ARRAY ARGUMENTS IN THE CALL STATEMENTS IN THIS ROUTINE.
159 ! IT HAS THE EFFECT OF PASSING ONLY THE FIRST ELEMENT OF THE ARRAY, RATHER
160 ! THAN THE ENTIRE ARRAY.  SEE:
161 ! http://www.mmm.ucar.edu/wrf/WG2/topics/deref_kludge.htm
163 !-----------------------------------------------------------------------
164 #include "deref_kludge.h"
165 !-----------------------------------------------------------------------
167 ! NEEDED BY SOME COMM LAYERS, E.G. RSL.  IF NEEDED, nmm_data_calls.inc IS
168 ! GENERATED FROM THE REGISTRY.  THE DEFINITION OF REGISTER_I1 ALLOWS
169 ! I1 DATA TO BE COMMUNICATED IN THIS ROUTINE IF NECESSARY.
171 !-----------------------------------------------------------------------
172 !-----------------------------------------------------------------------
173 !***********************************************************************
174 !***
175 !***               THE MAIN TIME INTEGRATION LOOP
176 !***
177 !-----------------------------------------------------------------------
179 !***  NTSD IS THE TIMESTEP COUNTER (Number of Time Steps Done)
181 !-----------------------------------------------------------------------
182 !***
183 !***  ADVANCE_count STARTS AT ZERO FOR ALL RUNS (REGULAR AND RESTART).
184 !***
185 !-----------------------------------------------------------------------
187       CALL DOMAIN_CLOCK_GET(GRID,ADVANCEcOUNT=NTSD_current)
189       IF(NTSD_current==0)THEN
190         IF(GRID%RESTART.AND.GRID%TSTART>0.)THEN
191           IHRST=NSTART_HOUR
192           NTSD_restart=NTSD
193         ELSE
194           IHRST=GRID%GMT
195           NSTART_HOUR=IHRST
196           NTSD_restart=0
197         ENDIF
198       ENDIF
200       NTSD=NTSD_restart+NTSD_current
201       LAST_TIME=domain_last_time_step(GRID)
203 !-----------------------------------------------------------------------
205 !!!!! IF(WRF_DM_ON_MONITOR() )THEN
206         WRITE(MESSAGE,125)NTSD,NTSD*GRID%DT/3600.
207   125   FORMAT(' SOLVE_NMM: TIMESTEP IS ',I5,'   TIME IS ',F7.3,' HOURS')
208         CALL WRF_MESSAGE(TRIM(MESSAGE))
209 !!!!  ENDIF
211 !-----------------------------------------------------------------------
212       CALL WRF_GET_DM_COMMUNICATOR(MPI_COMM_COMP)
213       CALL WRF_GET_NPROC(NPES)
214       CALL WRF_GET_MYPROC(MYPROC)
215       MYPE=MYPROC
216 !-----------------------------------------------------------------------
218 !***  OBTAIN DIMENSION INFORMATION STORED IN THE GRID DATA STRUCTURE.
220       CALL GET_IJK_FROM_GRID(GRID                                       &
221      &                      ,IDS,IDE,JDS,JDE,KDS,KDE                    &
222      &                      ,IMS,IME,JMS,JME,KMS,KME                    &
223      &                      ,IPS,IPE,JPS,JPE,KPS,KPE )
224 !-----------------------------------------------------------------------
226 !***  COMPUTE THESE STARTING AND STOPPING LOCATIONS FOR EACH TILE AND
227 !***  NUMBER OF TILES.
228 !***  SEE: http://www.mmm.ucar.edu/wrf/WG2/topics/settiles
230       CALL SET_TILES(GRID,IDS,IDE,JDS,JDE,IPS,IPE,JPS,JPE)
232 !-----------------------------------------------------------------------
233 !***  SET FLAG FOR THE OPERATIONAL PHYSICS SUITE.
234 !***  THIS WILL BE USED TO SAVE CLOCKTIME BY SKIPPING
235 !***  FREQUENT UPDATES OF THE MOIST ARRAY AND INSTEAD
236 !***  UPDATE IT ONLY WHEN IT IS NEEDED FOR PHYSICS.
237 !-----------------------------------------------------------------------
239       OPERATIONAL_PHYSICS=.FALSE.
241       IF(CONFIG_FLAGS%RA_SW_PHYSICS    ==GFDLSWSCHEME.AND.              &
242      &   CONFIG_FLAGS%RA_LW_PHYSICS    ==GFDLLWSCHEME.AND.              &
243      &   CONFIG_FLAGS%SF_SFCLAY_PHYSICS==MYJSFCSCHEME.AND.              &
244      &   CONFIG_FLAGS%BL_PBL_PHYSICS   ==MYJPBLSCHEME.AND.              &
245      &   CONFIG_FLAGS%CU_PHYSICS       ==BMJSCHEME.AND.                 &
246      &   CONFIG_FLAGS%MP_PHYSICS       ==ETAMPNEW)THEN
248         OPERATIONAL_PHYSICS=.TRUE.
250       ENDIF
252 !-----------------------------------------------------------------------
254 !***  TTEN, QTEN are used by GD convection scheme
256       ALLOCATE(TTEN(IMS:IME,KMS:KME,JMS:JME),STAT=ISTAT)
257       ALLOCATE(QTEN(IMS:IME,KMS:KME,JMS:JME),STAT=ISTAT)
258       ALLOCATE(RTHBLTEN(IMS:IME,KMS:KME,JMS:JME),STAT=ISTAT)
259       ALLOCATE(RQVBLTEN(IMS:IME,KMS:KME,JMS:JME),STAT=ISTAT)
260       ALLOCATE(RTHRATEN(IMS:IME,KMS:KME,JMS:JME),STAT=ISTAT)
261 #ifdef WRF_CHEM
262   ALLOCATE(CHEM_TRANS(IMS:IME,JMS:JME,KMS:KME,1:NUM_CHEM),STAT=ISTAT)
263   NUMGAS          = GET_LAST_GAS(CONFIG_FLAGS%CHEM_OPT)
264 #endif
266       IF(CONFIG_FLAGS%CU_PHYSICS==GDSCHEME)THEN
267         DO J=JMS,JME
268         DO K=KMS,KME
269         DO I=IMS,IME
270           TTEN(I,K,J)=T(I,J,K)
271           QTEN(I,K,J)=Q(I,J,K)
272         ENDDO
273         ENDDO
274         ENDDO
275       ENDIF
277       GRID%SIGMA=1 
278       HYDRO=.FALSE.
281       IDF=IDE-1
282       JDF=JDE-1
283       KDF=KDE-1
285 !-----------------------------------------------------------------------
287 !***  FOR NOW SET CONTROLS FOR TILES TO PATCHES
289 !-----------------------------------------------------------------------
290       ITS=IPS
291       ITE=MIN(IPE,IDF)
292       JTS=JPS
293       JTE=MIN(JPE,JDF)
294       KTS=KPS
295       KTE=MIN(KPE,KDF)
296       if(ntsd==0)then
297         write(0,*)' its=',its,' ite=',ite
298         write(0,*)' jts=',jts,' jte=',jte
299         write(0,*)' kts=',kts,' kte=',kte
300       endif
301 !-----------------------------------------------------------------------
302 !***  SET TIMING VARIABLES TO ZERO AT START OF FORECAST.
303 !-----------------------------------------------------------------------
304       if(ntsd==0)then
305         solve_tim=0.
306         exch_tim=0.
307         pdte_tim=0.
308         adve_tim=0.
309         vtoa_tim=0.
310         vadz_tim=0.
311         hadz_tim=0.
312         eps_tim=0.
313         vad2_tim=0.
314         had2_tim=0.
315         radiation_tim=0.
316         rdtemp_tim=0.
317         turbl_tim=0.
318         cltend_tim=0.
319         cucnvc_tim=0.
320         gsmdrive_tim=0.
321         hdiff_tim=0.
322         bocoh_tim=0.
323         pfdht_tim=0.
324         ddamp_tim=0.
325         bocov_tim=0.
326         uv_htov_tim=0.
327         exch_tim_max=0.
328         adjppt_tim=0.
329       endif
330 !-----------------------------------------------------------------------
331       N_MOIST=NUM_MOIST
333       DO J=MYJS_P4,MYJE_P4
334         IHEG(J)=MOD(J+1,2)
335         IHWG(J)=IHEG(J)-1
336         IVEG(J)=MOD(J,2)
337         IVWG(J)=IVEG(J)-1
338       ENDDO
340       DO J=MYJS_P4,MYJE_P4
341         IVW(J)=IVWG(J)
342         IVE(J)=IVEG(J)
343         IHE(J)=IHEG(J)
344         IHW(J)=IHWG(J)
345       ENDDO
347 !***  LATERAL POINTS IN THE BOUNDARY ARRAYS
349       LB=2*(IDF-IDS+1)+(JDF-JDS+1)-3
351 !***  APPROXIMATE GRIDPOINT SPACING (METERS)
353       JC=JMS+(JME-JMS)/2
354       GPS=SQRT(DX_NMM(IMS,JC)**2+DY_NMM**2)
356 !***  TIMESTEPS PER HOUR
358       TSPH=3600./GRID%DT
360       n_print_time=nint(3600./grid%dt)   ! Print stats once per hour
361 !-----------------------------------------------------------------------
363       NBOCO=0
365 !-----------------------------------------------------------------------
367 #if (NMM_NEST == 1)
368 !-----------------------------------------------------------------------------
369 !***  PATCHING NESTED BOUNDARIES.
370 !-----------------------------------------------------------------------------
372       CALL wrf_debug ( 100 , 'nmm: in patch' )
374       btimx=timef()
375 #ifdef DM_PARALLEL
376 #    include "HALO_NMM_ZZ.inc"
377 #endif
379       IF(GRID%ID/=1)THEN
381         CALL NESTBC_PATCH (PD_BXS,PD_BXE,PD_BYS,PD_BYE,T_BXS,T_BXE,T_BYS,T_BYE,Q_BXS,Q_BXE         &
382                           ,Q_BYS,Q_BYE,U_BXS,U_BXE,U_BYS,U_BYE,V_BXS,V_BXE,V_BYS,V_BYE             &
383                           ,Q2_BXS,Q2_BXE,Q2_BYS,Q2_BYE,CWM_BXS,CWM_BXE,CWM_BYS,CWM_BYE             &
384                           ,PD_BTXS,PD_BTXE,PD_BTYS,PD_BTYE,T_BTXS,T_BTXE,T_BTYS,T_BTYE             &
385                           ,Q_BTXS,Q_BTXE,Q_BTYS,Q_BTYE,U_BTXS,U_BTXE,U_BTYS,U_BTYE                 &
386                           ,V_BTXS,V_BTXE,V_BTYS,V_BTYE,Q2_BTXS,Q2_BTXE,Q2_BTYS,Q2_BTYE             &
387                           ,CWM_BTXS,CWM_BTXE,CWM_BTYS,CWM_BTYE,PDNEST_B,TNEST_B,QNEST_B,UNEST_B    &
388                           ,VNEST_B,Q2NEST_B,CWMNEST_B,PDNEST_BT,TNEST_BT,QNEST_BT                  &
389                           ,UNEST_BT,VNEST_BT,Q2NEST_BT,CWMNEST_BT                                  &
390                           ,GRID%SPEC_BDY_WIDTH                                                     &
391                           ,IDS,IDF,JDS,JDF,KDS,KDE                                                 &
392                           ,IMS,IME,JMS,JME,KMS,KME                                                 &
393                           ,ITS,ITE,JTS,JTE,KTS,KTE                                                 )
395         CALL wrf_debug ( 100 , 'nmm: out of patch' )
397 #ifdef MOVE_NESTS
399         IF(GRID%ID/=1.AND.MOD(NTSD,1)==0.AND.GRID%NUM_MOVES==-99)THEN
400           XLOC_1=(IDE-1)/2     ! This maneuvers the storm to the center of the nest quickly
401           YLOC_1=(JDE-1)/2     ! This maneuvers the storm to the center of the nest quickly
402         ENDIF
403 #endif
405       ENDIF
406 #endif
408 !-----------------------------------------------------------------------
409 !***  ALLOCATE PPTDAT ARRAY (PRECIP ASSIM):
410 !-----------------------------------------------------------------------
412       IF(GRID%PCPFLG.AND..NOT.ALLOCATED(PPTDAT))THEN
413         ALLOCATE(PPTDAT(IMS:IME,JMS:JME,3),STAT=ISTAT)
414       ENDIF
416 !-----------------------------------------------------------------------
417 !***
418 !***      Call READPCP to
419 !***            1) READ IN PRECIPITATION FOR HOURS 1, 2 and 3;
420 !***            2) Initialize DDATA to 999. (this is the amount
421 !***               of input precip allocated to each physics time step
422 !***               in ADJPPT; TURBL/SURFCE, which uses DDATA, is called
423 !***               before ADJPPT)
424 !***            3) Initialize LSPA to zero
425 !***
426 !-----------------------------------------------------------------------
427       IF (NTSD==0) THEN
428         IF (GRID%PCPFLG) THEN
429           CALL READPCP(PPTDAT,DDATA,LSPA                                &
430      &      ,IDS,IDE,JDS,JDE,KDS,KDE                                    &
431      &      ,IMS,IME,JMS,JME,KMS,KME                                    &
432      &      ,ITS,ITE,JTS,JTE,KTS,KTE)
433         ENDIF
434       ENDIF
435 !-----------------------------------------------------------------------
437       btim=timef()
439 !-----------------------------------------------------------------------
440 !***  ZERO OUT ACCUMULATED QUANTITIES WHEN NEEDED.
441 !-----------------------------------------------------------------------
443       CALL BUCKETS(NTSD,NPREC,NSRFC,NRDSW,NRDLW                         &
444      &            ,GRID%RESTART,GRID%TSTART                             &
445      &            ,NCLOD,NHEAT,GRID%NPHS,TSPH                           &
446      &            ,ACPREC,CUPREC,ACSNOW,ACSNOM,SSROFF,BGROFF            &
447      &            ,SFCEVP,POTEVP,SFCSHX,SFCLHX,SUBSHX,SNOPCX            &
448      &            ,SFCUVX,POTFLX                                        &
449      &            ,ARDSW,ASWIN,ASWOUT,ASWTOA                            &
450      &            ,ARDLW,ALWIN,ALWOUT,ALWTOA                            &
451      &            ,ACFRST,NCFRST,ACFRCV,NCFRCV                          &
452      &            ,AVCNVC,AVRAIN,TCUCN,TRAIN                            &
453      &            ,ASRFC                                                &
454      &            ,T,TLMAX,TLMIN,TSHLTR,PSHLTR,QSHLTR                   &
455      &            ,T02_MAX,T02_MIN,RH02_MAX,RH02_MIN                    &
456      &            ,IDS,IDE,JDS,JDE,KDS,KDE                              &
457      &            ,IMS,IME,JMS,JME,KMS,KME                              &
458      &            ,ITS,ITE,JTS,JTE,KTS,KTE)
459 !-----------------------------------------------------------------------
461       IF(NTSD==0)THEN
462         FIRST=.TRUE.
463 !       call hpm_init()
464         btimx=timef()
466 !-----------------------------------------------------------------------
467 #ifdef DM_PARALLEL
468 #    include "HALO_NMM_A.inc"
469 #endif
471 !-----------------------------------------------------------------------
472 #ifdef DM_PARALLEL
473       IF(CONFIG_FLAGS%MP_PHYSICS/=ETAMPNEW)THEN
474 #    include "HALO_NMM_A_3.inc"
475       ENDIF
476 #endif
477 !-----------------------------------------------------------------------
479 !***  Only for chemistry:
481 #ifdef WRF_CHEM   
482 #ifdef DM_PARALLEL
483 #    include "HALO_NMM_A_2.inc"
484 #endif
485 #endif
487 !-----------------------------------------------------------------------
488 !***  USE THE FOLLOWING VARIABLES TO KEEP TRACK OF EXCHANGE TIMES.
489 !-----------------------------------------------------------------------
490         exch_tim=exch_tim+timef()-btimx
491 !       this_tim=timef()-btimx
492 !       call mpi_allreduce(this_tim,et_max,1,mpi_real,mpi_max           &
493 !    &                    ,mpi_comm_comp,irtn)
494 !       exch_tim_max=exch_tim_max+et_max
495 !-----------------------------------------------------------------------
497         GO TO 2003
498       ENDIF
500 !-----------------------------------------------------------------------
501 !-----------------------------------------------------------------------
502  2000 CONTINUE
503 !-----------------------------------------------------------------------
504 !-----------------------------------------------------------------------
506 !-----------------------------------------------------------------------
507 !***  PRESSURE TENDENCY, SIGMA DOT, VERTICAL PART OF OMEGA-ALPHA
508 !-----------------------------------------------------------------------
510       btimx=timef()
511 !-----------------
512 #ifdef DM_PARALLEL
513 #    include "HALO_NMM_D.inc"
514 #endif
515 !-----------------
516       exch_tim=exch_tim+timef()-btimx
517 !     this_tim=timef()-btimx
518 !     call mpi_allreduce(this_tim,et_max,1,mpi_real,mpi_max             &
519 !    &                  ,mpi_comm_comp,irtn)
520 !     exch_tim_max=exch_tim_max+et_max
522       btimx=timef()
524       CALL PDTE(                                                        &
525 #ifdef DM_PARALLEL
526      &            GRID,MYPE,MPI_COMM_COMP,                              &
527 #endif
528      &            NTSD,GRID%DT,PT,ETA2,RES,HYDRO,HBM2                   &
529      &           ,PD,PDSL,PDSLO                                         &
530      &           ,PETDT,DIV,PSDT                                        &
531      &           ,IHE,IHW,IVE,IVW                                       &
532      &           ,IDS,IDF,JDS,JDF,KDS,KDE                               &
533      &           ,IMS,IME,JMS,JME,KMS,KME                               &
534      &           ,ITS,ITE,JTS,JTE,KTS,KTE)
536       pdte_tim=pdte_tim+timef()-btimx
538 !-----------------------------------------------------------------------
539 !***  ADVECTION OF T, U, AND V
540 !-----------------------------------------------------------------------
542       btimx=timef()
543 !-----------------
544 #ifdef DM_PARALLEL
545 #    include "HALO_NMM_F.inc"
546 #    include "HALO_NMM_F1.inc"
547 #endif
548 !-----------------
549       exch_tim=exch_tim+timef()-btimx
550 !     this_tim=timef()-btimx
551 !     call mpi_allreduce(this_tim,et_max,1,mpi_real,mpi_max             &
552 !    &                  ,mpi_comm_comp,irtn)
553 !     exch_tim_max=exch_tim_max+et_max
554       btimx=timef()
556       CALL ADVE(NTSD,GRID%DT,DETA1,DETA2,PDTOP                          &
557      &         ,CURV,F,FAD,F4D,EM_LOC,EMT_LOC,EN,ENT,DX_NMM,DY_NMM      &
558      &         ,HBM2,VBM2                                               &
559      &         ,T,U,V,PDSLO,TOLD,UOLD,VOLD                              &
560      &         ,PETDT,UPSTRM                                            &
561      &         ,FEW,FNS,FNE,FSE                                         &
562      &         ,ADT,ADU,ADV                                             & 
563      &         ,N_IUP_H,N_IUP_V                                         &
564      &         ,N_IUP_ADH,N_IUP_ADV                                     &
565      &         ,IUP_H,IUP_V,IUP_ADH,IUP_ADV                             &
566      &         ,IHE,IHW,IVE,IVW                                         &
567      &         ,IDS,IDF,JDS,JDF,KDS,KDE                                 &
568      &         ,IMS,IME,JMS,JME,KMS,KME                                 &
569      &         ,ITS,ITE,JTS,JTE,KTS,KTE)
571       adve_tim=adve_tim+timef()-btimx
573 !-----------------------------------------------------------------------
574 !***  PRESSURE TENDENCY, ETA/SIGMADOT, VERTICAL PART OF OMEGA-ALPHA TERM
575 !-----------------------------------------------------------------------
577       btimx=timef()
579       CALL VTOA(                                                        &
580 #ifdef DM_PARALLEL
581      &          GRID,                                                   &
582 #endif
583      &          NTSD,GRID%DT,PT,ETA2                                    &
584      &         ,HBM2,EF4T                                               &
585      &         ,T,DWDT,RTOP,OMGALF                                      &
586      &         ,PINT,DIV,PSDT,RES                                       &
587      &         ,IHE,IHW,IVE,IVW                                         &
588      &         ,IDS,IDF,JDS,JDF,KDS,KDE                                 &
589      &         ,IMS,IME,JMS,JME,KMS,KME                                 &
590      &         ,ITS,ITE,JTS,JTE,KTS,KTE)
592       vtoa_tim=vtoa_tim+timef()-btimx
594 !-----------------------------------------------------------------------
595 !***  VERTICAL ADVECTION OF HEIGHT
596 !-----------------------------------------------------------------------
598       btimx=timef()
600       CALL VADZ(NTSD,GRID%DT,FIS,GRID%SIGMA,DFL,HBM2                    &
601      &         ,DETA1,DETA2,PDTOP                                       &
602      &         ,PINT,PDSL,PDSLO,PETDT                                   &
603      &         ,RTOP,T,Q,CWM,Z,W,DWDT,PDWDT                             &
604      &         ,IHE,IHW,IVE,IVW                                         &
605      &         ,IDS,IDF,JDS,JDF,KDS,KDE                                 &
606      &         ,IMS,IME,JMS,JME,KMS,KME                                 &
607      &         ,ITS,ITE,JTS,JTE,KTS,KTE)
609       vadz_tim=vadz_tim+timef()-btimx
611 !-----------------------------------------------------------------------
612 !***  HORIZONTAL ADVECTION OF HEIGHT
613 !-----------------------------------------------------------------------
615       btimx=timef()
616 !-----------------
617 #ifdef DM_PARALLEL
618 #    include "HALO_NMM_G.inc"
619 #endif
620 !-----------------
621       exch_tim=exch_tim+timef()-btimx
622 !     this_tim=timef()-btimx
623 !     call mpi_allreduce(this_tim,et_max,1,mpi_real,mpi_max             &
624 !    &                  ,mpi_comm_comp,irtn)
625 !     exch_tim_max=exch_tim_max+et_max
627       btimx=timef()
629       CALL HADZ(NTSD,GRID%DT,HYDRO,HBM2,DETA1,DETA2,PDTOP               &
630      &         ,DX_NMM,DY_NMM,FAD                                       &
631      &         ,FEW,FNS,FNE,FSE                                         &
632      &         ,PDSL,U,V,W,Z                                            &
633      &         ,IHE,IHW,IVE,IVW                                         &
634      &         ,IDS,IDF,JDS,JDF,KDS,KDE                                 &
635      &         ,IMS,IME,JMS,JME,KMS,KME                                 &
636      &         ,ITS,ITE,JTS,JTE,KTS,KTE)
638       hadz_tim=hadz_tim+timef()-btimx
640 !-----------------------------------------------------------------------
641 !***  ADVECTION OF W
642 !-----------------------------------------------------------------------
644       btimx=timef()
645 !-----------------
646 #ifdef DM_PARALLEL
647 #    include "HALO_NMM_H.inc"
648 #endif
649 !-----------------
650       exch_tim=exch_tim+timef()-btimx
651 !     this_tim=timef()-btimx
652 !     call mpi_allreduce(this_tim,et_max,1,mpi_real,mpi_max             &
653 !    &                  ,mpi_comm_comp,irtn)
654 !     exch_tim_max=exch_tim_max+et_max
656       btimx=timef()
658       CALL EPS(NTSD,GRID%DT,HYDRO,DX_NMM,DY_NMM,FAD                     &
659      &        ,DETA1,DETA2,PDTOP,PT                                     &
660      &        ,HBM2,HBM3                                                &
661      &        ,PDSL,PDSLO,PINT,RTOP,PETDT,PDWDT                         &
662      &        ,DWDT,DWDTMN,DWDTMX                                       &
663      &        ,FNS,FEW,FNE,FSE                                          &
664      &        ,T,U,V,W,Q,CWM                                            &
665      &        ,IHE,IHW,IVE,IVW                                          &
666      &        ,IDS,IDF,JDS,JDF,KDS,KDE                                  &
667      &        ,IMS,IME,JMS,JME,KMS,KME                                  &
668      &        ,ITS,ITE,JTS,JTE,KTS,KTE)
670       eps_tim=eps_tim+timef()-btimx
672 !-----------------------------------------------------------------------
673 !***  VERTICAL ADVECTION OF Q, TKE, AND CLOUD WATER
674 !-----------------------------------------------------------------------
676       IF(MOD(NTSD,GRID%IDTAD)==0)THEN
677         btimx=timef()
679         vad2_micro_check: IF(CONFIG_FLAGS%MP_PHYSICS==ETAMPNEW)THEN
680           CALL VAD2(NTSD,GRID%DT,GRID%IDTAD,DX_NMM,DY_NMM               &
681      &             ,AETA1,AETA2,DETA1,DETA2,PDSL,PDTOP,HBM2             &
682      &             ,Q,Q2,CWM,PETDT                                      &
683      &             ,N_IUP_H,N_IUP_V                                     &
684      &             ,N_IUP_ADH,N_IUP_ADV                                 &
685      &             ,IUP_H,IUP_V,IUP_ADH,IUP_ADV                         &
686      &             ,IHE,IHW,IVE,IVW                                     &
687      &             ,IDS,IDF,JDS,JDF,KDS,KDE                             &
688      &             ,IMS,IME,JMS,JME,KMS,KME                             &
689      &             ,ITS,ITE,JTS,JTE,KTS,KTE)
691         ELSE vad2_micro_check
692           CALL VAD2_SCAL(NTSD,GRID%DT,GRID%IDTAD,DX_NMM,DY_NMM          &
693      &                  ,AETA1,AETA2,DETA1,DETA2,PDSL,PDTOP             &
694      &                  ,HBM2                                           &
695      &                  ,Q2,PETDT                                       &
696      &                  ,N_IUP_H,N_IUP_V                                &
697      &                  ,N_IUP_ADH,N_IUP_ADV                            &
698      &                  ,IUP_H,IUP_V,IUP_ADH,IUP_ADV                    &
699      &                  ,IHE,IHW,IVE,IVW                                &
700      &                  ,1,1                                            &
701      &                  ,IDS,IDF,JDS,JDF,KDS,KDE                        &
702      &                  ,IMS,IME,JMS,JME,KMS,KME                        &
703      &                  ,ITS,ITE,JTS,JTE,KTS,KTE)                              
704      
705           CALL VAD2_SCAL(NTSD,GRID%DT,GRID%IDTAD,DX_NMM,DY_NMM          &
706      &                  ,AETA1,AETA2,DETA1,DETA2,PDSL,PDTOP             &
707      &                  ,HBM2                                           &
708      &                  ,MOIST,PETDT                                    &
709      &                  ,N_IUP_H,N_IUP_V                                &
710      &                  ,N_IUP_ADH,N_IUP_ADV                            &
711      &                  ,IUP_H,IUP_V,IUP_ADH,IUP_ADV                    &
712      &                  ,IHE,IHW,IVE,IVW                                &
713      &                  ,NUM_MOIST,2                                    &
714      &                  ,IDS,IDF,JDS,JDF,KDS,KDE                        &
715      &                  ,IMS,IME,JMS,JME,KMS,KME                        &
716      &                  ,ITS,ITE,JTS,JTE,KTS,KTE)
718           CALL VAD2_SCAL(NTSD,GRID%DT,GRID%IDTAD,DX_NMM,DY_NMM          &
719      &                  ,AETA1,AETA2,DETA1,DETA2,PDSL,PDTOP             &
720      &                  ,HBM2                                           &
721      &                  ,SCALAR,PETDT                                   &
722      &                  ,N_IUP_H,N_IUP_V                                &
723      &                  ,N_IUP_ADH,N_IUP_ADV                            &
724      &                  ,IUP_H,IUP_V,IUP_ADH,IUP_ADV                    &
725      &                  ,IHE,IHW,IVE,IVW                                &
726      &                  ,NUM_SCALAR,2                                   &
727      &                  ,IDS,IDF,JDS,JDF,KDS,KDE                        &
728      &                  ,IMS,IME,JMS,JME,KMS,KME                        &
729      &                  ,ITS,ITE,JTS,JTE,KTS,KTE)
732           DO K=KTS,KTE
733           DO J=MYJS,MYJE
734           DO I=MYIS,MYIE
735             Q(I,J,K)=MOIST(I,J,K,P_QV)/(1.+MOIST(I,J,K,P_QV))
736           ENDDO
737           ENDDO   
738           ENDDO   
740         ENDIF vad2_micro_check
742         vad2_tim=vad2_tim+timef()-btimx
744       ENDIF
745 !    
746 !-----------------------------------------------------------------------
747 !***  VERTICAL ADVECTION OF CHEMISTRY
748 !-----------------------------------------------------------------------
750 #ifdef WRF_CHEM  
751       IF(MOD(NTSD,GRID%IDTAD)==0)THEN
752 #ifdef IBM
753         btimx=timef()
754 #endif
756         DO L=1,NUM_CHEM
757         DO K=KMS,KME
758         DO J=JMS,JME
759         DO I=IMS,IME
760            CHEM_TRANS(I,J,K,L)=CHEM(I,K,J,L)
761         ENDDO
762         ENDDO
763         ENDDO
764         ENDDO
766          
767         CALL VAD2_SCAL(NTSD,GRID%DT,GRID%IDTAD,DX_NMM,DY_NMM            &
768      &           ,AETA1,AETA2,DETA1,DETA2,PDSL,PDTOP                    &
769      &           ,HBM2                                                  &
770      &           ,CHEM_TRANS,PETDT                                            & 
771      &           ,N_IUP_H,N_IUP_V                                       &
772      &           ,N_IUP_ADH,N_IUP_ADV                                   &
773      &           ,IUP_H,IUP_V,IUP_ADH,IUP_ADV                           &
774      &           ,IHE,IHW,IVE,IVW                                       &
775      &           ,NUM_CHEM,1                                            & 
776      &           ,IDS,IDF,JDS,JDF,KDS,KDE                               &
777      &           ,IMS,IME,JMS,JME,KMS,KME                               &
778      &           ,ITS,ITE,JTS,JTE,KTS,KTE)
780       ENDIF      
781 #endif 
783 !-----------------------------------------------------------------------
784 !***  HORIZONTAL ADVECTION OF Q, TKE, AND CLOUD WATER
785 !-----------------------------------------------------------------------
787       IF(MOD(NTSD,GRID%IDTAD)==0)THEN
788         btimx=timef()
789 !-----------------
790 #ifdef DM_PARALLEL
791 #    include "HALO_NMM_I.inc"
792 #endif
794 #ifdef DM_PARALLEL
795         IF(CONFIG_FLAGS%MP_PHYSICS/=ETAMPNEW)THEN
796 #    include "HALO_NMM_I_3.inc"
797         ENDIF
798 #endif
800 !-----------------
801         exch_tim=exch_tim+timef()-btimx
802 !       this_tim=timef()-btimx
803 !       call mpi_allreduce(this_tim,et_max,1,mpi_real,mpi_max           &
804 !    &                    ,mpi_comm_comp,irtn)
805 !       exch_tim_max=exch_tim_max+et_max
807         btimx=timef()
809 !-----------------------------------------------------------------------
810         had2_micro_check: IF(CONFIG_FLAGS%MP_PHYSICS==ETAMPNEW)THEN
811 !-----------------------------------------------------------------------
813           CALL HAD2(                                                   &
814 #if defined(DM_PARALLEL)
815      &              GRID%DOMDESC,                                      &
816 #endif
817      &              NTSD,GRID%DT,GRID%IDTAD,DX_NMM,DY_NMM              &
818      &             ,AETA1,AETA2,DETA1,DETA2,PDSL,PDTOP                 &
819      &             ,HBM2,HBM3                                          &
820      &             ,Q,Q2,CWM,U,V,Z,HYDRO                               &
821      &             ,N_IUP_H,N_IUP_V                                    &
822      &             ,N_IUP_ADH,N_IUP_ADV                                &
823      &             ,IUP_H,IUP_V,IUP_ADH,IUP_ADV                        &
824      &             ,IHE,IHW,IVE,IVW                                    &
825      &             ,IDS,IDF,JDS,JDF,KDS,KDE                            &
826      &             ,IMS,IME,JMS,JME,KMS,KME                            &
827      &             ,ITS,ITE,JTS,JTE,KTS,KTE)
829 !***  UPDATE MOIST ARRAY.
830 !***  REMEMBER THAT MOIST IS ONLY USED WITH THE PHYSICS AND THUS
831 !***  THE P_QV SLOT (=2) IS MIXING RATIO, NOT SPECIFIC HUMIDITY.
832 !***  ALTHOUGH MOIST IS ONLY USED FOR PHYSICS IN OPERATIONS, IT IS 
833 !***  UPDATED HERE FROM Q EVERY ADVECTION TIMESTEP FOR NON-OPERATIONAL
834 !***  CONFIGURATIONS WHERE IT MAY BE USED OUTSIDE OF THE PHYSICS.
836           IF(.NOT.OPERATIONAL_PHYSICS)THEN
837             DO K=KTS,KTE
838             DO J=MYJS,MYJE
839             DO I=MYIS,MYIE
840               MOIST(I,J,K,P_QV)=Q(I,J,K)/(1.-Q(I,J,K))
841               WC = CWM(I,J,K)
842               QI = 0.
843               QR = 0.
844               QW = 0.
845               FICE=F_ICE(I,K,J)
846               FRAIN=F_RAIN(I,K,J)
848               IF(FICE>=1.)THEN
849                 QI=WC
850               ELSEIF(FICE<=0.)THEN
851                 QW=WC
852               ELSE
853                 QI=FICE*WC
854                 QW=WC-QI
855               ENDIF
857               IF(QW>0..AND.FRAIN>0.)THEN
858                 IF(FRAIN>=1.)THEN
859                   QR=QW
860                   QW=0.
861                 ELSE
862                   QR=FRAIN*QW
863                   QW=QW-QR
864                 ENDIF
865               ENDIF
867               MOIST(I,J,K,P_QC)=QW
868               MOIST(I,J,K,P_QR)=QR
869               MOIST(I,J,K,P_QI)=0.
870               MOIST(I,J,K,P_QS)=QI
871               MOIST(I,J,K,P_QG)=0.
872             ENDDO
873             ENDDO
874             ENDDO
875           ENDIF
877 !-----------------------------------------------------------------------
878         ELSE had2_micro_check
879 !-----------------------------------------------------------------------
881           CALL HAD2_SCAL(                                              &
882 #if defined(DM_PARALLEL)
883      &                   GRID%DOMDESC,                                 &
884 #endif        
885      &                   NTSD,GRID%DT,GRID%IDTAD,DX_NMM,DY_NMM         &
886      &                  ,AETA1,AETA2,DETA1,DETA2,PDSL,PDTOP            &
887      &                  ,HBM2,HBM3                                     &
888      &                  ,Q2,U,V,Z,HYDRO                                &
889      &                  ,N_IUP_H,N_IUP_V                               &
890      &                  ,N_IUP_ADH,N_IUP_ADV                           &
891      &                  ,IUP_H,IUP_V,IUP_ADH,IUP_ADV                   &
892      &                  ,IHE,IHW,IVE,IVW                               &
893      &                  ,1,1                                           &
894      &                  ,IDS,IDF,JDS,JDF,KDS,KDE                       &
895      &                  ,IMS,IME,JMS,JME,KMS,KME                       &
896      &                  ,ITS,ITE,JTS,JTE,KTS,KTE)
897 !        
898           CALL HAD2_SCAL(                                              &
899 #if defined(DM_PARALLEL)
900      &                   GRID%DOMDESC,                                 &
901 #endif
902      &                   NTSD,GRID%DT,GRID%IDTAD,DX_NMM,DY_NMM         &
903      &                  ,AETA1,AETA2,DETA1,DETA2,PDSL,PDTOP            &
904      &                  ,HBM2,HBM3                                     &
905      &                  ,MOIST,U,V,Z,HYDRO                             &
906      &                  ,N_IUP_H,N_IUP_V                               &
907      &                  ,N_IUP_ADH,N_IUP_ADV                           &
908      &                  ,IUP_H,IUP_V,IUP_ADH,IUP_ADV                   &
909      &                  ,IHE,IHW,IVE,IVW                               &
910      &                  ,NUM_MOIST,2                                   &
911      &                  ,IDS,IDF,JDS,JDF,KDS,KDE                       &
912      &                  ,IMS,IME,JMS,JME,KMS,KME                       &
913      &                  ,ITS,ITE,JTS,JTE,KTS,KTE)
915           CALL HAD2_SCAL(                                              &
916 #if defined(DM_PARALLEL)
917      &                   GRID%DOMDESC,                                 &
918 #endif
919      &                   NTSD,GRID%DT,GRID%IDTAD,DX_NMM,DY_NMM         &
920      &                  ,AETA1,AETA2,DETA1,DETA2,PDSL,PDTOP            &
921      &                  ,HBM2,HBM3                                     &
922      &                  ,SCALAR,U,V,Z,HYDRO                            &
923      &                  ,N_IUP_H,N_IUP_V                               &
924      &                  ,N_IUP_ADH,N_IUP_ADV                           &
925      &                  ,IUP_H,IUP_V,IUP_ADH,IUP_ADV                   &
926      &                  ,IHE,IHW,IVE,IVW                               &
927      &                  ,NUM_SCALAR,2                                  &
928      &                  ,IDS,IDF,JDS,JDF,KDS,KDE                       &
929      &                  ,IMS,IME,JMS,JME,KMS,KME                       &
930      &                  ,ITS,ITE,JTS,JTE,KTS,KTE)                             
931 !    
932           DO K=KTS,KTE 
933           DO J=MYJS,MYJE
934           DO I=MYIS,MYIE
935             Q(I,J,K)=MOIST(I,J,K,P_QV)/(1.+MOIST(I,J,K,P_QV))           
936           ENDDO
937           ENDDO    
938           ENDDO   
940 !-----------------------------------------------------------------------
941         ENDIF had2_micro_check
942 !-----------------------------------------------------------------------
944         had2_tim=had2_tim+timef()-btimx
945       ENDIF
947 !-----------------------------------------------------------------------
948 !***  HORIZONTAL ADVECTION OF CHEMISTRY
949 !-----------------------------------------------------------------------
951 #ifdef WRF_CHEM
952       IF(MOD(NTSD,GRID%IDTAD)==0)THEN
953         btimx=timef()
954 #ifdef DM_PARALLEL
955 #    include "HALO_NMM_I_2.inc"
956 #endif
957         exch_tim=exch_tim+timef()-btimx
958 !       this_tim=timef()-btimx
959 !       call mpi_allreduce(this_tim,et_max,1,mpi_real,mpi_max           &
960 !    &                    ,mpi_comm_comp,irtn)
961 !       exch_tim_max=exch_tim_max+et_max
963         btimx=timef()
965         CALL HAD2_SCAL(                                                &
966 #if defined(DM_PARALLEL)
967      &                 GRID%DOMDESC,                                   &
968 #endif
969      &                 NTSD,GRID%DT,GRID%IDTAD,DX_NMM,DY_NMM           &
970      &                ,AETA1,AETA2,DETA1,DETA2,PDSL,PDTOP              &
971      &                ,HBM2,HBM3                                       &
972      &                ,CHEM_TRANS,U,V,Z,HYDRO                                &
973      &                ,N_IUP_H,N_IUP_V                                 &
974      &                ,N_IUP_ADH,N_IUP_ADV                             &
975      &                ,IUP_H,IUP_V,IUP_ADH,IUP_ADV                     &
976      &                ,IHE,IHW,IVE,IVW                                 &
977      &                ,NUM_CHEM,1                                      &
978      &                ,IDS,IDF,JDS,JDF,KDS,KDE                         &
979      &                ,IMS,IME,JMS,JME,KMS,KME                         &
980      &                ,ITS,ITE,JTS,JTE,KTS,KTE)
982       ENDIF
983         DO L=1,NUM_CHEM
984         DO J=JMS,JME
985         DO K=KMS,KME
986         DO I=IMS,IME
987            CHEM(I,K,J,L)=CHEM_TRANS(I,J,K,L)
988         ENDDO
989         ENDDO
990         ENDDO
991         ENDDO
992 #endif
995 !----------------------------------------------------------------------
996 !***  RADIATION
997 !----------------------------------------------------------------------
999 !***  When allocating CAM radiation 4d arrays (ozmixm, aerosolc), 
1000 !***  the following two scalars are not needed.
1002       NUM_OZMIXM=1
1003       NUM_AEROSOLC=1
1005       IF(NTSD<=0)THEN
1006         NTSD_rad=NTSD
1007       ELSE
1009 !***  Call radiation just BEFORE the top of the hour
1010 !***  so that updated fields are written to history files.
1012         NTSD_rad=NTSD+1
1013       ENDIF
1015       IF(MOD(NTSD_rad,GRID%NRADS)==0.OR.                               &
1016      &   MOD(NTSD_rad,GRID%NRADL)==0)THEN
1018         btimx=timef()
1019         IF(OPERATIONAL_PHYSICS)THEN
1020           CALL UPDATE_MOIST(MOIST,Q,CWM,F_ICE,F_RAIN,N_MOIST           &
1021      &                     ,IDS,IDF,JDS,JDF,KDS,KDE                    &
1022      &                     ,IMS,IME,JMS,JME,KMS,KME                    &
1023      &                     ,ITS,ITE,JTS,JTE,KTS,KTE)
1024         ENDIF
1026         CALL RADIATION(NTSD_rad,GRID%DT,GRID%JULDAY,GRID%JULYR         &
1027      &                ,GRID%XTIME,GRID%JULIAN                          &
1028      &                ,IHRST,GRID%NPHS                                 &
1029      &                ,GLAT,GLON,GRID%NRADS,GRID%NRADL                 &
1030      &                ,DETA1,DETA2,AETA1,AETA2,ETA1,ETA2,PDTOP,PT      &
1031      &                ,PD,RES,PINT,T,Q,MOIST,THS,ALBEDO,EPSR           &
1032      &                ,F_ICE,F_RAIN                                    &
1033 #ifdef WRF_CHEM
1034      &                ,GD_CLOUD,GD_CLOUD2                              &
1035 #endif
1036      &                ,SM,HBM2,CLDFRA,N_MOIST,RESTRT                   &
1037      &                ,RLWTT,RSWTT,RLWIN,RSWIN,RSWINC,RSWOUT           &
1038      &                ,RLWTOA,RSWTOA,CZMEAN                            &
1039      &                ,CFRACL,CFRACM,CFRACH,SIGT4                      &
1040      &                ,ACFRST,NCFRST,ACFRCV,NCFRCV                     &
1041      &                ,CUPPT,VEGFRC,SNO,HTOP,HBOT                      &
1042      &                ,Z,SICE,NUM_AEROSOLC,NUM_OZMIXM                  &
1043      &                ,GRID,CONFIG_FLAGS                               &
1044      &                ,RTHRATEN                                        &  
1045 #ifdef WRF_CHEM
1046      &                ,PM2_5_DRY, PM2_5_WATER, PM2_5_DRY_EC            &
1047      &                ,TAUAER1, TAUAER2, TAUAER3, TAUAER4              & 
1048      &                ,GAER1, GAER2, GAER3, GAER4                      &
1049      &                ,WAER1, WAER2, WAER3, WAER4                      &
1050 #endif
1051      &                ,IDS,IDF,JDS,JDF,KDS,KDE                         &
1052      &                ,IMS,IME,JMS,JME,KMS,KME                         &
1053      &                ,ITS,ITE,JTS,JTE,KTS,KTE)
1055         DO J=JMS,JME
1056         DO I=IMS,IME
1057           GSW(I,J)=RSWIN(I,J)-RSWOUT(I,J)
1058         ENDDO
1059         ENDDO
1061 !                            ***  NOTE  ***
1062 ! RLWIN/RSWIN  - downward longwave/shortwave at the surface (formerly TOTLWDN/TOTSWDN)
1063 ! RSWINC - CLEAR-SKY downward shortwave at the surface (new for AQ)
1064 !                            ***  NOTE  ***
1066         radiation_tim=radiation_tim+timef()-btimx
1067       ENDIF
1069 !----------------------------------------------------------------------
1070 !***  APPLY TEMPERATURE TENDENCY DUE TO RADIATION
1071 !----------------------------------------------------------------------
1073       btimx=timef()
1075 !     Pass in XTIME (elapsed time from start of parent) to compute
1076 !     the time passed into the zenith angle code consistently between
1077 !     RDTEMP and RADIATION.
1079       CALL RDTEMP(NTSD,GRID%DT,GRID%JULDAY,GRID%JULYR                  &
1080      &           ,GRID%XTIME,IHRST,GLAT,GLON                           &
1081      &           ,CZEN,CZMEAN,T,RSWTT,RLWTT,HBM2                       &
1082      &           ,IDS,IDF,JDS,JDF,KDS,KDE                              &
1083      &           ,IMS,IME,JMS,JME,KMS,KME                              &
1084      &           ,ITS,ITE,JTS,JTE,KTS,KTE)
1086       rdtemp_tim=rdtemp_tim+timef()-btimx
1088 !----------------------------------------------------------------------
1089 !***  TURBULENT PROCESSES 
1090 !----------------------------------------------------------------------
1092       IF(MOD(NTSD,GRID%NPHS)==0)THEN
1094         btimx=timef()
1096         IF(OPERATIONAL_PHYSICS                                         &
1097      &    .AND.MOD(NTSD_rad,GRID%NRADS)/=0                             &
1098      &    .AND.MOD(NTSD_rad,GRID%NRADL)/=0)THEN
1099           CALL UPDATE_MOIST(MOIST,Q,CWM,F_ICE,F_RAIN,N_MOIST           &
1100      &                     ,IDS,IDF,JDS,JDF,KDS,KDE                    &
1101      &                     ,IMS,IME,JMS,JME,KMS,KME                    &
1102      &                     ,ITS,ITE,JTS,JTE,KTS,KTE)
1103         ENDIF
1105         CALL TURBL(NTSD,GRID%DT,GRID%NPHS,RESTRT                       &
1106      &            ,N_MOIST,GRID%NUM_SOIL_LAYERS,SLDPTH,DZSOIL          &
1107      &            ,DETA1,DETA2,AETA1,AETA2,ETA1,ETA2,PDTOP,PT          &
1108      &            ,SM,HBM2,VBM2,DX_NMM,DFRLG                           &
1109      &            ,CZEN,CZMEAN,SIGT4,RLWIN,RSWIN,RADOT                 &
1110      &            ,PD,RES,PINT,T,Q,CWM,F_ICE,F_RAIN,SR                 &
1111      &            ,Q2,U,V,THS,NMM_TSK,SST,PREC,SNO                     &
1112      &            ,FIS,Z0,Z0BASE,USTAR,PBLH,LPBL,EL_MYJ                &
1113      &            ,MOIST,RMOL,MOL                                      &
1114      &            ,EXCH_H,AKHS,AKMS,AKHS_OUT,AKMS_OUT                  &
1115      &            ,THZ0,QZ0,UZ0,VZ0,QSH,MAVAIL                         &
1116      &            ,STC,SMC,CMC,SMSTAV,SMSTOT,SSROFF,BGROFF             &
1117      &            ,IVGTYP,ISLTYP,VEGFRC,SHDMIN,SHDMAX,GRNFLX           &
1118      &            ,SFCEXC,ACSNOW,ACSNOM,SNOPCX,SICE,TG,SOILTB          &
1119      &            ,ALBASE,MXSNAL,ALBEDO,SH2O,SI,EPSR,EMBCK             &
1120      &            ,U10,V10,TH10,Q10,TSHLTR,QSHLTR,PSHLTR               &
1121      &            ,T2,QSG,QVG,QCG,SOILT1,TSNAV,SMFR3D,KEEPFR3DFLAG     &
1122      &            ,TWBS,QWBS,SFCSHX,SFCLHX,SFCEVP                      &
1123      &            ,POTEVP,POTFLX,SUBSHX                                &
1124      &            ,APHTIM,ARDSW,ARDLW,ASRFC                            &
1125      &            ,RSWOUT,RSWTOA,RLWTOA                                &
1126      &            ,ASWIN,ASWOUT,ASWTOA,ALWIN,ALWOUT,ALWTOA             &
1127      &            ,UZ0H,VZ0H,DUDT,DVDT                                 & 
1128      &            ,RTHBLTEN,RQVBLTEN                                   & 
1129      &            ,GRID%PCPFLG,DDATA                                   &
1130      &            ,GRID,CONFIG_FLAGS                                   &
1131      &            ,IHE,IHW,IVE,IVW                                     &
1132      &            ,IDS,IDF,JDS,JDF,KDS,KDE                             &
1133      &            ,IMS,IME,JMS,JME,KMS,KME                             &
1134      &            ,ITS,ITE,JTS,JTE,KTS,KTE)
1136 !                     ***  NOTE  ***
1137 ! RLWIN/RSWIN - downward longwave/shortwave at the surface
1138 !                     ***  NOTE  ***
1140         turbl_tim=turbl_tim+timef()-btimx
1142         btimx=timef()
1143 !-----------------
1144 #ifdef DM_PARALLEL
1145 # include "HALO_NMM_TURBL_A.inc"
1146 #endif
1147 !-----------------
1148 #ifdef DM_PARALLEL
1149 # include "HALO_NMM_TURBL_B.inc"
1150 #endif
1151 !-----------------
1152         exch_tim=exch_tim+timef()-btimx
1153 !       this_tim=timef()-btimx
1154 !       call mpi_allreduce(this_tim,et_max,1,mpi_real,mpi_max           &
1155 !    &                    ,mpi_comm_comp,irtn)
1156 !       exch_tim_max=exch_tim_max+et_max
1158 !***  INTERPOLATE WINDS FROM H POINTS BACK TO V POINTS.
1160         btimx=timef()
1161         CALL UV_H_TO_V(NTSD,GRID%DT,GRID%NPHS,UZ0H,VZ0H,UZ0,VZ0         &
1162      &                ,DUDT,DVDT,U,V,HBM2,IVE,IVW                       &
1163      &                ,IDS,IDF,JDS,JDF,KDS,KDE                          &
1164      &                ,IMS,IME,JMS,JME,KMS,KME                          &
1165      &                ,ITS,ITE,JTS,JTE,KTS,KTE)
1166         uv_htov_tim=uv_htov_tim+timef()-btimx
1168 !----------------------------------------------------------------------
1169 !*** STORE ORIGINAL TEMPERATURE ARRAY
1170 !----------------------------------------------------------------------
1172         btimx=timef()
1173 !-----------------
1174 #ifdef DM_PARALLEL
1175 #    include "HALO_NMM_J.inc"
1176 #endif
1178 #ifdef DM_PARALLEL
1179         IF(CONFIG_FLAGS%MP_PHYSICS/=ETAMPNEW)THEN
1180 #    include "HALO_NMM_J_3.inc"
1181         ENDIF
1182 #endif
1184 #ifdef WRF_CHEM
1185 #ifdef DM_PARALLEL
1186 #    include "HALO_NMM_J_2.inc"
1187 #endif
1188 #endif
1189 !-----------------
1190         exch_tim=exch_tim+timef()-btimx
1191 !       this_tim=timef()-btimx
1192 !       call mpi_allreduce(this_tim,et_max,1,mpi_real,mpi_max           &
1193 !    &                    ,mpi_comm_comp,irtn)
1194 !       exch_tim_max=exch_tim_max+et_max
1196         ICLTEND=-1
1197         btimx=timef()
1199         CALL CLTEND(ICLTEND,GRID%NPHS,T,T_OLD,T_ADJ                    &
1200      &             ,IDS,IDF,JDS,JDF,KDS,KDE                            &
1201      &             ,IMS,IME,JMS,JME,KMS,KME                            &
1202      &             ,ITS,ITE,JTS,JTE,KTS,KTE)
1204         cltend_tim=cltend_tim+timef()-btimx
1205       ENDIF
1207 !----------------------------------------------------------------------
1208 !***  CONVECTIVE PRECIPITATION
1209 !----------------------------------------------------------------------
1211       IF(MOD(NTSD,GRID%NCNVC)==0.AND.                                  &
1212      &   CONFIG_FLAGS%CU_PHYSICS==KFETASCHEME)THEN
1214         btimx=timef()
1215 !-----------------
1216 #ifdef DM_PARALLEL
1217 #    include "HALO_NMM_C.inc"
1218 #endif
1219 !-----------------
1220         exch_tim=exch_tim+timef()-btimx
1221 !       this_tim=timef()-btimx
1222 !       call mpi_allreduce(this_tim,et_max,1,mpi_real,mpi_max          &
1223 !    &                    ,mpi_comm_comp,irtn)
1224 !       exch_tim_max=exch_tim_max+et_max
1225       ENDIF
1227       convection: IF(CONFIG_FLAGS%CU_PHYSICS/=0)THEN
1228         btimx=timef()
1230 !***  GET TENDENCIES FOR GD SCHEME.
1231 !    
1232         IF(CONFIG_FLAGS%CU_PHYSICS==GDSCHEME)THEN
1233           DT_INV=1./GRID%DT
1234           DO J=JMS,JME
1235           DO K=KMS,KME
1236           DO I=IMS,IME
1237             TTEN(I,K,J)=(T(I,J,K)-TTEN(I,K,J))*DT_INV
1238             QTEN(I,K,J)=(Q(I,J,K)-QTEN(I,K,J))*DT_INV
1239           ENDDO
1240           ENDDO
1241           ENDDO
1242         ENDIF
1244 !***  UPDATE THE MOIST ARRAY IF NEEDED.
1246         IF(OPERATIONAL_PHYSICS                                         &
1247      &    .AND.MOD(NTSD_rad,GRID%NRADS)/=0                             &
1248      &    .AND.MOD(NTSD_rad,GRID%NRADL)/=0                             &
1249      &    .AND.MOD(NTSD,GRID%NPHS)/=0)THEN
1250           CALL UPDATE_MOIST(MOIST,Q,CWM,F_ICE,F_RAIN,N_MOIST           &
1251      &                     ,IDS,IDF,JDS,JDF,KDS,KDE                    &
1252      &                     ,IMS,IME,JMS,JME,KMS,KME                    &
1253      &                     ,ITS,ITE,JTS,JTE,KTS,KTE)
1254         ENDIF
1256 !----------------------------------------------------------------------
1257         CALL CUCNVC(NTSD,GRID%DT,GRID%NCNVC,GRID%NRADS,GRID%NRADL      &
1258      &             ,GPS,RESTRT,HYDRO,CLDEFI,N_MOIST,GRID%ENSDIM        &
1259      &             ,MOIST                                              &
1260      &             ,DETA1,DETA2,AETA1,AETA2,ETA1,ETA2                  &
1261      &             ,F_ICE,F_RAIN                                       &
1262 !***  Changes for other cu schemes, most for GD scheme
1263      &             ,APR_GR,APR_W,APR_MC,TTEN,QTEN                      &
1264      &             ,APR_ST,APR_AS,APR_CAPMA                            &
1265      &             ,APR_CAPME,APR_CAPMI                                &
1266      &             ,MASS_FLUX,XF_ENS                                   &
1267      &             ,PR_ENS,GSW                                         &
1268 #ifdef WRF_CHEM
1269      &             ,GD_CLOUD,GD_CLOUD2,RAINCV                          &
1270 #endif
1272      &             ,PDTOP,PT,PD,RES,PINT,T,Q,CWM,TCUCN                 &
1273      &             ,OMGALF,U,V,W,Z,FIS,W0AVG                           &
1274      &             ,PREC,ACPREC,CUPREC,CUPPT,CPRATE                    &
1275      &             ,SM,HBM2,LPBL,CNVBOT,CNVTOP                         &
1276      &             ,HTOP,HBOT,HTOPD,HBOTD,HTOPS,HBOTS                  &
1277      &             ,RTHBLTEN,RQVBLTEN,RTHRATEN                         & 
1278      &             ,AVCNVC,ACUTIM,IHE,IHW                              &
1279      &             ,GRID,CONFIG_FLAGS                                  &
1280      &             ,IDS,IDF,JDS,JDF,KDS,KDE                            &
1281      &             ,IMS,IME,JMS,JME,KMS,KME                            &
1282      &             ,IPS,IPE,JPS,JPE,KPS,KPE                            &
1283      &             ,ITS,ITE,JTS,JTE,KTS,KTE)
1284 !----------------------------------------------------------------------
1286         cucnvc_tim=cucnvc_tim+timef()-btimx
1288       ENDIF convection
1290 !----------------------------------------------------------------------
1291 !***  GRIDSCALE MICROPHYSICS (CONDENSATION & PRECIPITATION)
1292 !----------------------------------------------------------------------
1294       IF(MOD(NTSD,GRID%NPHS)==0)THEN
1295         btimx=timef()
1297         CALL GSMDRIVE(NTSD,GRID%DT,GRID%NPHS,N_MOIST                   &
1298      &               ,DX_NMM(ITS,JC),GRID%DY,SM,HBM2,FIS               &
1299      &               ,DETA1,DETA2,AETA1,AETA2,ETA1,ETA2                &
1300      &               ,PDTOP,PT,PD,RES,PINT,T,Q,CWM,TRAIN               &
1301      &               ,MOIST,SCALAR,NUM_SCALAR                          &
1302      &               ,F_ICE,F_RAIN,F_RIMEF,SR                          &
1303      &               ,PREC,ACPREC,AVRAIN                               &
1304      &               ,MP_RESTART_STATE                                 &
1305      &               ,TBPVS_STATE                                      &
1306      &               ,TBPVS0_STATE                                     &
1307      &               ,GRID,CONFIG_FLAGS                                &
1308      &               ,IDS,IDF,JDS,JDF,KDS,KDE                          &
1309      &               ,IMS,IME,JMS,JME,KMS,KME                          &
1310      &               ,ITS,ITE,JTS,JTE,KTS,KTE)
1312         gsmdrive_tim=gsmdrive_tim+timef()-btimx
1314 !-----------------------------------------------------------------------
1315 !---------PRECIPITATION ASSIMILATION------------------------------------
1316 !-----------------------------------------------------------------------
1318         IF (GRID%PCPFLG) THEN
1319           btimx=timef()
1321           CALL CHKSNOW(NTSD,GRID%DT,GRID%NPHS,SR,PPTDAT                 &
1322      &      ,IDS,IDE,JDS,JDE,KDS,KDE                                    &
1323      &      ,IMS,IME,JMS,JME,KMS,KME                                    &
1324      &      ,ITS,ITE,JTS,JTE,KTS,KTE)
1325           CALL ADJPPT(NTSD,GRID%DT,GRID%NPHS,PREC,LSPA,PPTDAT,DDATA     &
1326      &      ,IDS,IDE,JDS,JDE,KDS,KDE                                    &
1327      &      ,IMS,IME,JMS,JME,KMS,KME                                    &
1328      &      ,ITS,ITE,JTS,JTE,KTS,KTE)
1330           adjppt_tim=adjppt_tim+timef()-btimx
1331         ENDIF
1333 !----------------------------------------------------------------------
1334 !***  CALCULATE TEMP TENDENCIES AND RESTORE ORIGINAL TEMPS
1335 !----------------------------------------------------------------------
1336 !      
1337         ICLTEND=0
1338         btimx=timef()
1340         CALL CLTEND(ICLTEND,GRID%NPHS,T,T_OLD,T_ADJ                    &
1341      &             ,IDS,IDF,JDS,JDF,KDS,KDE                            &
1342      &             ,IMS,IME,JMS,JME,KMS,KME                            &
1343      &             ,ITS,ITE,JTS,JTE,KTS,KTE)
1345         cltend_tim=cltend_tim+timef()-btimx
1346       ENDIF
1348 !----------------------------------------------------------------------
1349 !***  UPDATE TEMP TENDENCIES FROM CLOUD PROCESSES EVERY TIME STEP
1350 !----------------------------------------------------------------------
1352       ICLTEND=1
1353       btimx=timef()
1355       CALL CLTEND(ICLTEND,GRID%NPHS,T,T_OLD,T_ADJ                      &
1356      &           ,IDS,IDF,JDS,JDF,KDS,KDE                              &
1357      &           ,IMS,IME,JMS,JME,KMS,KME                              &
1358      &           ,ITS,ITE,JTS,JTE,KTS,KTE)
1360       cltend_tim=cltend_tim+timef()-btimx
1362 !----------------------------------------------------------------------
1363 !***  LATERAL DIFFUSION
1364 !----------------------------------------------------------------------
1366       btimx=timef()
1367 !-----------------
1368 #ifdef DM_PARALLEL
1369 #    include "HALO_NMM_K.inc"
1370 #endif
1371 !-----------------
1372       exch_tim=exch_tim+timef()-btimx
1373 !     this_tim=timef()-btimx
1374 !     call mpi_allreduce(this_tim,et_max,1,mpi_real,mpi_max             &
1375 !    &                  ,mpi_comm_comp,irtn)
1376 !     exch_tim_max=exch_tim_max+et_max
1378       btimx=timef()
1380       CALL HDIFF(NTSD,GRID%DT,FIS,DY_NMM,HDAC,HDACV                    &
1381      &          ,HBM2,DETA1,GRID%SIGMA                                 &
1382      &          ,T,Q,U,V,Q2,Z,W,SM,SICE                                &
1383      &          ,IHE,IHW,IVE,IVW                                       &
1384      &          ,IDS,IDF,JDS,JDF,KDS,KDE                               &
1385      &          ,IMS,IME,JMS,JME,KMS,KME                               &
1386      &          ,ITS,ITE,JTS,JTE,KTS,KTE)
1388       IF(.NOT.OPERATIONAL_PHYSICS)THEN
1389         DO K=KTS,KTE
1390         DO J=MYJS,MYJE
1391         DO I=MYIS,MYIE
1392 !!!       MOIST(I,J,K,P_QV)=MAX(0.,Q(I,J,K)/(1.-Q(I,J,K)))
1393           MOIST(I,J,K,P_QV)=Q(I,J,K)/(1.-Q(I,J,K))           !<-- Update mixing ratio
1394         ENDDO
1395         ENDDO
1396         ENDDO
1397       ENDIF
1399       hdiff_tim=hdiff_tim+timef()-btimx
1401 !----------------------------------------------------------------------
1402 !***  UPDATING BOUNDARY VALUES AT HEIGHT POINTS
1403 !----------------------------------------------------------------------
1405       btimx=timef()
1406 !-----------------
1407 #ifdef DM_PARALLEL
1408 #    include "HALO_NMM_L.inc"
1409 #endif
1411 #ifdef DM_PARALLEL
1412 #    include "HALO_NMM_L_3.inc"
1413 #endif
1415 #ifdef WRF_CHEM
1416 #ifdef DM_PARALLEL
1417 #    include "HALO_NMM_L_2.inc"
1418 #endif
1419 #endif
1420 !-----------------
1421       exch_tim=exch_tim+timef()-btimx
1422 !     this_tim=timef()-btimx
1423 !     call mpi_allreduce(this_tim,et_max,1,mpi_real,mpi_max             &
1424 !    &                  ,mpi_comm_comp,irtn)
1425 !     exch_tim_max=exch_tim_max+et_max
1427       btimx=timef()
1429       CALL BOCOH(GRID%ID,NTSD,GRID%DT,NEST,NUNIT_NBC,NBOCO,LAST_TIME,TSPH &
1430      &          ,LB,ETA1,ETA2,PDTOP,PT,RES                              &
1431      &          ,PD_BXS,PD_BXE,PD_BYS,PD_BYE,T_BXS,T_BXE,T_BYS,T_BYE    &
1432      &          ,Q_BXS,Q_BXE,Q_BYS,Q_BYE,U_BXS,U_BXE,U_BYS,U_BYE,V_BXS  &
1433      &          ,V_BXE,V_BYS,V_BYE,Q2_BXS,Q2_BXE,Q2_BYS,Q2_BYE,CWM_BXS  &
1434      &          ,CWM_BXE,CWM_BYS,CWM_BYE,PD_BTXS,PD_BTXE,PD_BTYS        &
1435      &          ,PD_BTYE,T_BTXS,T_BTXE,T_BTYS,T_BTYE,Q_BTXS,Q_BTXE      &
1436      &          ,Q_BTYS,Q_BTYE,U_BTXS,U_BTXE,U_BTYS,U_BTYE,V_BTXS       &
1437      &          ,V_BTXE,V_BTYS,V_BTYE,Q2_BTXS,Q2_BTXE,Q2_BTYS,Q2_BTYE   &
1438      &          ,CWM_BTXS,CWM_BTXE,CWM_BTYS,CWM_BTYE,PD,T,Q,Q2,CWM,PINT &
1439      &          ,MOIST,N_MOIST,SCALAR,NUM_SCALAR                        &
1440 #ifdef WRF_CHEM
1441      &          ,CHEM,NUMGAS,CONFIG_FLAGS                               &
1442 #endif
1443      &          ,GRID%SPEC_BDY_WIDTH,Z                                  &
1444      &          ,IHE,IHW,IVE,IVW                                        &
1445      &          ,IDS,IDF,JDS,JDF,KDS,KDE                                &
1446      &          ,IMS,IME,JMS,JME,KMS,KME                                &
1447      &          ,ITS,ITE,JTS,JTE,KTS,KTE)
1451       bocoh_tim=bocoh_tim+timef()-btimx
1452 !     if(mod(ntsd,n_print_time)==0)then
1453 !       call twr(t,0,'t',ntsd,mype,npes,mpi_comm_comp &
1454 !    &          ,ids,ide,jds,jde,kds,kde                               &
1455 !    &          ,ims,ime,jms,jme,kms,kme                               &
1456 !    &          ,its,ite,jts,jte,kts,kte)
1457 !     endif
1459 !----------------------------------------------------------------------
1460 !***  IS IT TIME FOR A CHECK POINT ON THE MODEL HISTORY FILE?
1461 !----------------------------------------------------------------------
1463  2003 CONTINUE
1465 !----------------------------------------------------------------------
1466 !***  PRESSURE GRD, CORIOLIS, DIVERGENCE, AND HORIZ PART OF OMEGA-ALPHA
1467 !----------------------------------------------------------------------
1469       btimx=timef()
1470 !-----------------
1471 #ifdef DM_PARALLEL
1472 #    include "HALO_NMM_A.inc"
1473 #endif
1475 #ifdef DM_PARALLEL
1476       IF(CONFIG_FLAGS%MP_PHYSICS/=ETAMPNEW)THEN
1477 #    include "HALO_NMM_A_3.inc"
1478       ENDIF
1479 #endif
1481 #ifdef WRF_CHEM
1482 #ifdef DM_PARALLEL
1483 #    include "HALO_NMM_A_2.inc"
1484 #endif
1485 #endif
1486 !-----------------
1487       exch_tim=exch_tim+timef()-btimx
1488 !     this_tim=timef()-btimx
1489 !     call mpi_allreduce(this_tim,et_max,1,mpi_real,mpi_max             &
1490 !    &                  ,mpi_comm_comp,irtn)
1491 !     exch_tim_max=exch_tim_max+et_max
1493       btimx=timef()
1495       CALL PFDHT(NTSD,LAST_TIME,PT,DETA1,DETA2,PDTOP,RES,FIS           &
1496      &          ,HYDRO,GRID%SIGMA,FIRST,DX_NMM,DY_NMM                  &
1497      &          ,HBM2,VBM2,VBM3                                        &
1498      &          ,FDIV,FCP,WPDAR,DFL,CPGFU,CPGFV                        &
1499      &          ,PD,PDSL,T,Q,U,V,CWM,OMGALF,PINT,DWDT                  &
1500      &          ,RTOP,DIV,FEW,FNS,FNE,FSE                              &
1501      &          ,IHE,IHW,IVE,IVW                                       &
1502      &          ,IDS,IDF,JDS,JDF,KDS,KDE                               &
1503      &          ,IMS,IME,JMS,JME,KMS,KME                               &
1504      &          ,ITS,ITE,JTS,JTE,KTS,KTE)
1507       pfdht_tim=pfdht_tim+timef()-btimx
1509 !----------------------------------------------------------------------
1510 !***  DIVERGENCE DAMPING
1511 !----------------------------------------------------------------------
1513       btimx=timef()
1514 !-----------------
1515 #ifdef DM_PARALLEL
1516 #    include "HALO_NMM_B.inc"
1517 #endif
1518 !-----------------
1519       exch_tim=exch_tim+timef()-btimx
1520 !     this_tim=timef()-btimx
1521 !     call mpi_allreduce(this_tim,et_max,1,mpi_real,mpi_max            &
1522 !    &                  ,mpi_comm_comp,irtn)
1523 !     exch_tim_max=exch_tim_max+et_max
1525       btimx=timef()
1527       CALL DDAMP(NTSD,GRID%DT,DETA1,DETA2,PDSL,PDTOP,DIV,HBM2          &
1528      &          ,T,U,V,DDMPU,DDMPV                                     &
1529      &          ,IHE,IHW,IVE,IVW                                       &
1530      &          ,IDS,IDF,JDS,JDF,KDS,KDE                               &
1531      &          ,IMS,IME,JMS,JME,KMS,KME                               &
1532      &          ,ITS,ITE,JTS,JTE,KTS,KTE)
1534       ddamp_tim=ddamp_tim+timef()-btimx
1536 !----------------------------------------------------------------------
1537 !----------------------------------------------------------------------
1539       IF(FIRST.AND.NTSD==0)THEN
1540         FIRST=.FALSE.
1541         btimx=timef()
1542 !-----------------
1543 #ifdef DM_PARALLEL
1544 #    include "HALO_NMM_A.inc"
1545 #endif
1546 #ifdef WRF_CHEM
1547 #ifdef DM_PARALLEL
1548 #    include "HALO_NMM_A_2.inc"
1549 #endif
1550 #endif
1551 !-----------------
1552         exch_tim=exch_tim+timef()-btimx
1553 !       this_tim=timef()-btimx
1554 !       call mpi_allreduce(this_tim,et_max,1,mpi_real,mpi_max          &
1555 !    &                    ,mpi_comm_comp,irtn)
1556 !       exch_tim_max=exch_tim_max+et_max
1558         GO TO 2000
1559       ENDIF
1561 !----------------------------------------------------------------------
1562 !***  UPDATING BOUNDARY VALUES AT VELOCITY POINTS
1563 !----------------------------------------------------------------------
1565       btimx=timef()
1566 !-----------------
1567 #ifdef DM_PARALLEL
1568 #    include "HALO_NMM_C.inc"
1569 #endif
1570 !-----------------
1571       exch_tim=exch_tim+timef()-btimx
1572 !     this_tim=timef()-btimx
1573 !     call mpi_allreduce(this_tim,et_max,1,mpi_real,mpi_max            &
1574 !    &                  ,mpi_comm_comp,irtn)
1575 !     exch_tim_max=exch_tim_max+et_max
1577       btimx=timef()
1579       CALL BOCOV(GRID%ID,NTSD,GRID%DT,LB,U_BXS,U_BXE,U_BYS,U_BYE,V_BXS &
1580      &          ,V_BXE,V_BYS,V_BYE,U_BTXS,U_BTXE,U_BTYS,U_BTYE,V_BTXS  &
1581      &          ,V_BTXE,V_BTYS,V_BTYE,U,V                              &
1582      &          ,GRID%SPEC_BDY_WIDTH                                   &
1583      &          ,IHE,IHW,IVE,IVW                                       &
1584      &          ,IDS,IDF,JDS,JDF,KDS,KDE                               &
1585      &          ,IMS,IME,JMS,JME,KMS,KME                               &
1586      &          ,ITS,ITE,JTS,JTE,KTS,KTE )
1589       bocov_tim=bocov_tim+timef()-btimx
1591 !----------------------------------------------------------------------
1592 !***  COPY THE NMM VARIABLE Q2 TO THE WRF VARIABLE TKE_MYJ
1593 !----------------------------------------------------------------------
1595       DO K=KTS,KTE
1596       DO J=JTS,JTE
1597       DO I=ITS,ITE
1598         TKE_MYJ(I,J,K)=0.5*Q2(I,J,K) !TKE is q squared over 2
1599       ENDDO
1600       ENDDO
1601       ENDDO
1603 !----------------------------------------------------------------------
1605       IF(LAST_TIME.AND.ALLOCATED(PPTDAT))THEN
1606         DEALLOCATE(PPTDAT,STAT=ISTAT)
1607       ENDIF
1609 !----------------------------------------------------------------------
1611       solve_tim=solve_tim+timef()-btim
1613 !----------------------------------------------------------------------
1614 !***  PRINT TIMING VARIABLES WHEN DESIRED.
1615 !----------------------------------------------------------------------
1617       sum_tim=pdte_tim+adve_tim+vtoa_tim+vadz_tim+hadz_tim+eps_tim     &
1618      &       +vad2_tim+had2_tim+radiation_tim+rdtemp_tim+turbl_tim     &
1619      &       +cltend_tim+cucnvc_tim+gsmdrive_tim+hdiff_tim             &
1620      &       +bocoh_tim+pfdht_tim+ddamp_tim+bocov_tim+uv_htov_tim      &
1621      &       +exch_tim+adjppt_tim
1623       if(mod(ntsd,n_print_time)==0)then
1624         write(0,*)' ntsd=',ntsd,' solve_tim=',solve_tim*1.e-3          &
1625      &           ,' sum_tim=',sum_tim*1.e-3
1626         write(0,*)' pdte_tim=',pdte_tim*1.e-3,' pct=',pdte_tim/sum_tim*100.
1627         write(0,*)' adve_tim=',adve_tim*1.e-3,' pct=',adve_tim/sum_tim*100.
1628         write(0,*)' vtoa_tim=',vtoa_tim*1.e-3,' pct=',vtoa_tim/sum_tim*100.
1629         write(0,*)' vadz_tim=',vadz_tim*1.e-3,' pct=',vadz_tim/sum_tim*100.
1630         write(0,*)' hadz_tim=',hadz_tim*1.e-3,' pct=',hadz_tim/sum_tim*100.
1631         write(0,*)' eps_tim=',eps_tim*1.e-3,' pct=',eps_tim/sum_tim*100.
1632         write(0,*)' vad2_tim=',vad2_tim*1.e-3,' pct=',vad2_tim/sum_tim*100.
1633         write(0,*)' had2_tim=',had2_tim*1.e-3,' pct=',had2_tim/sum_tim*100.
1634         write(0,*)' radiation_tim=',radiation_tim*1.e-3,' pct=',radiation_tim/sum_tim*100.
1635         write(0,*)' rdtemp_tim=',rdtemp_tim*1.e-3,' pct=',rdtemp_tim/sum_tim*100.
1636         write(0,*)' turbl_tim=',turbl_tim*1.e-3,' pct=',turbl_tim/sum_tim*100.
1637         write(0,*)' cltend_tim=',cltend_tim*1.e-3,' pct=',cltend_tim/sum_tim*100.
1638         write(0,*)' cucnvc_tim=',cucnvc_tim*1.e-3,' pct=',cucnvc_tim/sum_tim*100.
1639         write(0,*)' gsmdrive_tim=',gsmdrive_tim*1.e-3,' pct=',gsmdrive_tim/sum_tim*100.
1640         write(0,*)' adjppt_tim=',adjppt_tim*1.e-3,' pct=',adjppt_tim/sum_tim*100.
1641         write(0,*)' hdiff_tim=',hdiff_tim*1.e-3,' pct=',hdiff_tim/sum_tim*100.
1642         write(0,*)' bocoh_tim=',bocoh_tim*1.e-3,' pct=',bocoh_tim/sum_tim*100.
1643         write(0,*)' pfdht_tim=',pfdht_tim*1.e-3,' pct=',pfdht_tim/sum_tim*100.
1644         write(0,*)' ddamp_tim=',ddamp_tim*1.e-3,' pct=',ddamp_tim/sum_tim*100.
1645         write(0,*)' bocov_tim=',bocov_tim*1.e-3,' pct=',bocov_tim/sum_tim*100.
1646         write(0,*)' uv_h_to_v_tim=',uv_htov_tim*1.e-3,' pct=',uv_htov_tim/sum_tim*100.
1647         write(0,*)' exch_tim=',exch_tim*1.e-3,' pct=',exch_tim/sum_tim*100.
1648 !        call time_stats(exch_tim,'exchange',ntsd,mype,npes,mpi_comm_comp)
1649 !        write(0,*)' exch_tim_max=',exch_tim_max*1.e-3
1651         call field_stats(t,mype,mpi_comm_comp                          &
1652      &                  ,ids,ide,jds,jde,kds,kde                       &
1653      &                  ,ims,ime,jms,jme,kms,kme                       &
1654      &                  ,its,ite,jts,jte,kts,kte)
1655       endif
1657 !     if(last_time)then
1658       DEALLOCATE(TTEN,STAT=ISTAT)
1659       DEALLOCATE(QTEN,STAT=ISTAT)
1660       DEALLOCATE(RTHRATEN,STAT=ISTAT)
1661       DEALLOCATE(RTHBLTEN,STAT=ISTAT)
1662       DEALLOCATE(RQVBLTEN,STAT=ISTAT)
1663 #ifdef WRF_CHEM
1664       DEALLOCATE(CHEM_TRANS,STAT=ISTAT)
1665 #endif
1667 #define COPY_OUT
1668 #include <scalar_derefs.inc>
1669       Return
1670 !----------------------------------------------------------------------
1671 !**********************************************************************
1672 !**********************************************************************
1673 !*************    EXIT FROM THE TIME LOOP    **************************
1674 !**********************************************************************
1675 !**********************************************************************
1676 !----------------------------------------------------------------------
1677       END SUBROUTINE SOLVE_NMM
1678 !----------------------------------------------------------------------
1679 !**********************************************************************
1680 !----------------------------------------------------------------------
1681       SUBROUTINE TWR(ARRAY,KK,FIELD,NTSD,MYPE,NPES,MPI_COMM_COMP       &
1682      &              ,IDS,IDE,JDS,JDE,KDS,KDE                           &
1683      &              ,IMS,IME,JMS,JME,KMS,KME                           &
1684      &              ,ITS,ITE,JTS,JTE,KTS,KTE)
1685 !----------------------------------------------------------------------
1686 !**********************************************************************
1687       USE MODULE_EXT_INTERNAL
1689       IMPLICIT NONE
1690 #if defined(DM_PARALLEL) && !defined(STUBMPI)
1691       INCLUDE "mpif.h"
1692 #endif
1693 !----------------------------------------------------------------------
1694 !----------------------------------------------------------------------
1695       INTEGER,INTENT(IN) :: IDS,IDE,JDS,JDE,KDS,KDE                    &
1696      &                     ,IMS,IME,JMS,JME,KMS,KME                    &
1697      &                     ,ITS,ITE,JTS,JTE,KTS,KTE                    &
1698      &                     ,KK,MPI_COMM_COMP,MYPE,NPES,NTSD
1700       REAL,DIMENSION(IMS:IME,JMS:JME,KMS:KME+KK),INTENT(IN) :: ARRAY
1702       CHARACTER(*),INTENT(IN) :: FIELD
1704 !*** LOCAL VARIABLES
1706 #if defined(DM_PARALLEL) && !defined(STUBMPI)
1707       INTEGER,DIMENSION(MPI_STATUS_SIZE) :: JSTAT
1708       INTEGER,DIMENSION(MPI_STATUS_SIZE,4) :: STATUS_ARRAY
1709 #endif
1710       INTEGER,DIMENSION(2) :: IM_REM,JM_REM,IT_REM,JT_REM
1712       INTEGER :: I,IENDX,IER,IPE,IRECV,IRTN,ISEND,IUNIT                &
1713      &          ,J,K,N,NLEN,NSIZE
1714       INTEGER :: ITS_REM,ITE_REM,JTS_REM,JTE_REM
1716       REAL,DIMENSION(IDS:IDE,JDS:JDE) :: TWRITE
1717       REAL,ALLOCATABLE,DIMENSION(:) :: VALUES
1718       CHARACTER(5) :: TIMESTEP
1719       CHARACTER(6) :: FMT
1720       CHARACTER(12) :: FILENAME
1721 !----------------------------------------------------------------------
1722 !**********************************************************************
1723 !----------------------------------------------------------------------
1725       IF(NTSD<=9)THEN
1726         FMT='(I1.1)'
1727         NLEN=1
1728       ELSEIF(NTSD<=99)THEN
1729         FMT='(I2.2)'
1730         NLEN=2
1731       ELSEIF(NTSD<=999)THEN
1732         FMT='(I3.3)'
1733         NLEN=3
1734       ELSEIF(NTSD<=9999)THEN
1735         FMT='(I4.4)'
1736         NLEN=4
1737       ELSEIF(NTSD<=99999)THEN
1738         FMT='(I5.5)'
1739         NLEN=5
1740       ENDIF
1741       WRITE(TIMESTEP,FMT)NTSD
1742       FILENAME=FIELD//'_'//TIMESTEP(1:NLEN)
1744       IF(MYPE==0)THEN
1745         CALL INT_GET_FRESH_HANDLE(IUNIT)
1746         CLOSE(IUNIT)
1747         OPEN(UNIT=IUNIT,FILE=FILENAME,FORM='UNFORMATTED',IOSTAT=IER)
1748       ENDIF
1750 !----------------------------------------------------------------------
1751 !!!!  DO 500 K=KTS,KTE+KK     !Unflipped
1752 !!!!  DO 500 K=KTE+KK,KTS,-1
1753       DO 500 K=KDE-1,KDS,-1   !Write LM layers top down for checking
1754 !----------------------------------------------------------------------
1756 #if defined(DM_PARALLEL) && !defined(STUBMPI)
1757       IF(MYPE==0)THEN
1758         DO J=JTS,JTE
1759         DO I=ITS,ITE
1760           TWRITE(I,J)=ARRAY(I,J,K)
1761         ENDDO
1762         ENDDO
1764         DO IPE=1,NPES-1
1765           CALL MPI_RECV(IT_REM,2,MPI_INTEGER,IPE,IPE                    &
1766      &                 ,MPI_COMM_COMP,JSTAT,IRECV)
1767           CALL MPI_RECV(JT_REM,2,MPI_INTEGER,IPE,IPE                    &
1768      &                 ,MPI_COMM_COMP,JSTAT,IRECV)
1770           ITS_REM=IT_REM(1)
1771           ITE_REM=IT_REM(2)
1772           JTS_REM=JT_REM(1)
1773           JTE_REM=JT_REM(2)
1775           NSIZE=(ITE_REM-ITS_REM+1)*(JTE_REM-JTS_REM+1)
1776           ALLOCATE(VALUES(1:NSIZE))
1778           CALL MPI_RECV(VALUES,NSIZE,MPI_REAL,IPE,IPE                   &
1779      &                 ,MPI_COMM_COMP,JSTAT,IRECV)
1780           N=0
1781           DO J=JTS_REM,JTE_REM
1782             DO I=ITS_REM,ITE_REM
1783               N=N+1
1784               TWRITE(I,J)=VALUES(N)
1785             ENDDO
1786           ENDDO
1788           DEALLOCATE(VALUES)
1790         ENDDO
1792 !----------------------------------------------------------------------
1793       ELSE
1794         NSIZE=(ITE-ITS+1)*(JTE-JTS+1)
1795         ALLOCATE(VALUES(1:NSIZE))
1797         N=0
1798         DO J=JTS,JTE
1799         DO I=ITS,ITE
1800           N=N+1
1801           VALUES(N)=ARRAY(I,J,K)
1802         ENDDO
1803         ENDDO
1805         IT_REM(1)=ITS
1806         IT_REM(2)=ITE
1807         JT_REM(1)=JTS
1808         JT_REM(2)=JTE
1810         CALL MPI_SEND(IT_REM,2,MPI_INTEGER,0,MYPE                       &
1811      &               ,MPI_COMM_COMP,ISEND)
1812         CALL MPI_SEND(JT_REM,2,MPI_INTEGER,0,MYPE                       &
1813      &               ,MPI_COMM_COMP,ISEND)
1815         CALL MPI_SEND(VALUES,NSIZE,MPI_REAL,0,MYPE                      &
1816      &               ,MPI_COMM_COMP,ISEND)
1818         DEALLOCATE(VALUES)
1820       ENDIF
1821 !----------------------------------------------------------------------
1823       CALL MPI_BARRIER(MPI_COMM_COMP,IRTN)
1825       IF(MYPE==0)THEN
1827         DO J=JDS,JDE-1
1828           IENDX=IDE-1
1829           IF(MOD(J,2)==0)IENDX=IENDX-1
1830           WRITE(IUNIT)(TWRITE(I,J),I=1,IENDX)
1831         ENDDO
1833       ENDIF
1834 #else
1836         DO J=JDS,JDE-1
1837           IENDX=IDE-1
1838           IF(MOD(J,2)==0)IENDX=IENDX-1
1839           WRITE(IUNIT)(ARRAY(I,K,J),I=1,IENDX)
1840         ENDDO
1842 #endif
1845 !----------------------------------------------------------------------
1846   500 CONTINUE
1848       IF(MYPE==0)CLOSE(IUNIT)
1849 !----------------------------------------------------------------------
1851       END SUBROUTINE TWR
1852 !----------------------------------------------------------------------
1853 !**********************************************************************
1854 !----------------------------------------------------------------------
1855       SUBROUTINE VWR(ARRAY,KK,FIELD,NTSD,MYPE,NPES,MPI_COMM_COMP       &
1856      &              ,IDS,IDE,JDS,JDE,KDS,KDE                           &
1857      &              ,IMS,IME,JMS,JME,KMS,KME                           &
1858      &              ,ITS,ITE,JTS,JTE,KTS,KTE)
1859 !----------------------------------------------------------------------
1860 !**********************************************************************
1861       USE MODULE_EXT_INTERNAL
1863       IMPLICIT NONE
1864 #if defined(DM_PARALLEL) && !defined(STUBMPI)
1865       INCLUDE "mpif.h"
1866 #endif
1867 !----------------------------------------------------------------------
1868 !----------------------------------------------------------------------
1869       INTEGER,INTENT(IN) :: IDS,IDE,JDS,JDE,KDS,KDE                    &
1870      &                     ,IMS,IME,JMS,JME,KMS,KME                    &
1871      &                     ,ITS,ITE,JTS,JTE,KTS,KTE                    &
1872      &                     ,KK,MPI_COMM_COMP,MYPE,NPES,NTSD
1874       REAL,DIMENSION(IMS:IME,JMS:JME,KMS:KME+KK),INTENT(IN) :: ARRAY
1876       CHARACTER(*),INTENT(IN) :: FIELD
1878 !*** LOCAL VARIABLES
1880 #if defined(DM_PARALLEL) && !defined(STUBMPI)
1881       INTEGER,DIMENSION(MPI_STATUS_SIZE) :: JSTAT
1882       INTEGER,DIMENSION(MPI_STATUS_SIZE,4) :: STATUS_ARRAY
1883 #endif
1884       INTEGER,DIMENSION(2) :: IM_REM,JM_REM,IT_REM,JT_REM
1886       INTEGER :: I,IENDX,IER,IPE,IRECV,IRTN,ISEND,IUNIT                &
1887      &          ,J,K,L,N,NLEN,NSIZE
1888       INTEGER :: ITS_REM,ITE_REM,JTS_REM,JTE_REM
1890       REAL,DIMENSION(IDS:IDE,JDS:JDE) :: TWRITE
1891       REAL,ALLOCATABLE,DIMENSION(:) :: VALUES
1892       CHARACTER(5) :: TIMESTEP
1893       CHARACTER(6) :: FMT
1894       CHARACTER(12) :: FILENAME
1895       LOGICAL :: OPENED
1896 !----------------------------------------------------------------------
1897 !**********************************************************************
1898 !----------------------------------------------------------------------
1900       IF(NTSD<=9)THEN
1901         FMT='(I1.1)'
1902         NLEN=1
1903       ELSEIF(NTSD<=99)THEN
1904         FMT='(I2.2)'
1905         NLEN=2
1906       ELSEIF(NTSD<=999)THEN
1907         FMT='(I3.3)'
1908         NLEN=3
1909       ELSEIF(NTSD<=9999)THEN
1910         FMT='(I4.4)'
1911         NLEN=4
1912       ELSEIF(NTSD<=99999)THEN
1913         FMT='(I5.5)'
1914         NLEN=5
1915       ENDIF
1916       WRITE(TIMESTEP,FMT)NTSD
1917       FILENAME=FIELD//'_'//TIMESTEP(1:NLEN)
1919       IF(MYPE==0)THEN
1920         CALL INT_GET_FRESH_HANDLE(IUNIT)
1921         CLOSE(IUNIT)
1922         OPEN(UNIT=IUNIT,FILE=FILENAME,FORM='UNFORMATTED',IOSTAT=IER)
1923       ENDIF
1924 !     IF(MYPE==0)THEN
1925 !       OPEN_UNIT: DO L=51,99
1926 !         INQUIRE(L,OPENED=OPENED)
1927 !         IF(.NOT.OPENED)THEN
1928 !           IUNIT=L
1929 !           OPEN(UNIT=IUNIT,FILE=FILENAME,STATUS='NEW'                 &
1930 !               ,FORM='UNFORMATTED',IOSTAT=IER)
1931 !           IF(IER/=0)THEN
1932 !             WRITE(0,*)' Opening file error=',IER
1933 !             WRITE(6,*)' Opening file error=',IER
1934 !           ENDIF
1935 !           EXIT OPEN_UNIT
1936 !         ENDIF
1937 !       ENDDO OPEN_UNIT
1938 !     ENDIF
1940 !----------------------------------------------------------------------
1941 !!!!  DO 500 K=KTS,KTE+KK     !Unflipped
1942 !!!!  DO 500 K=KTE+KK,KTS,-1
1943       DO 500 K=KDE-1,KDS,-1   !Write LM layers top down for checking
1944 !----------------------------------------------------------------------
1946 #if defined(DM_PARALLEL) && !defined(STUBMPI)
1947       IF(MYPE==0)THEN
1948         DO J=JTS,JTE
1949         DO I=ITS,ITE
1950           TWRITE(I,J)=ARRAY(I,J,K)
1951         ENDDO
1952         ENDDO
1954         DO IPE=1,NPES-1
1955           CALL MPI_RECV(IT_REM,2,MPI_INTEGER,IPE,IPE                    &
1956      &                 ,MPI_COMM_COMP,JSTAT,IRECV)
1957           CALL MPI_RECV(JT_REM,2,MPI_INTEGER,IPE,IPE                    &
1958      &                 ,MPI_COMM_COMP,JSTAT,IRECV)
1960           ITS_REM=IT_REM(1)
1961           ITE_REM=IT_REM(2)
1962           JTS_REM=JT_REM(1)
1963           JTE_REM=JT_REM(2)
1965           NSIZE=(ITE_REM-ITS_REM+1)*(JTE_REM-JTS_REM+1)
1966           ALLOCATE(VALUES(1:NSIZE))
1968           CALL MPI_RECV(VALUES,NSIZE,MPI_REAL,IPE,IPE                   &
1969      &                 ,MPI_COMM_COMP,JSTAT,IRECV)
1970           N=0
1971           DO J=JTS_REM,JTE_REM
1972             DO I=ITS_REM,ITE_REM
1973               N=N+1
1974               TWRITE(I,J)=VALUES(N)
1975             ENDDO
1976           ENDDO
1978           DEALLOCATE(VALUES)
1980         ENDDO
1982 !----------------------------------------------------------------------
1983       ELSE
1984         NSIZE=(ITE-ITS+1)*(JTE-JTS+1)
1985         ALLOCATE(VALUES(1:NSIZE))
1987         N=0
1988         DO J=JTS,JTE
1989         DO I=ITS,ITE
1990           N=N+1
1991           VALUES(N)=ARRAY(I,J,K)
1992         ENDDO
1993         ENDDO
1995         IT_REM(1)=ITS
1996         IT_REM(2)=ITE
1997         JT_REM(1)=JTS
1998         JT_REM(2)=JTE
2000         CALL MPI_SEND(IT_REM,2,MPI_INTEGER,0,MYPE                       &
2001      &               ,MPI_COMM_COMP,ISEND)
2002         CALL MPI_SEND(JT_REM,2,MPI_INTEGER,0,MYPE                       &
2003      &               ,MPI_COMM_COMP,ISEND)
2005         CALL MPI_SEND(VALUES,NSIZE,MPI_REAL,0,MYPE                      &
2006      &               ,MPI_COMM_COMP,ISEND)
2008         DEALLOCATE(VALUES)
2010       ENDIF
2011 !----------------------------------------------------------------------
2013       CALL MPI_BARRIER(MPI_COMM_COMP,IRTN)
2015       IF(MYPE==0)THEN
2017         DO J=JDS,JDE-1
2018           IENDX=IDE-1
2019           IF(MOD(J,2)==1)IENDX=IENDX-1
2020           WRITE(IUNIT)(TWRITE(I,J),I=1,IENDX)
2021         ENDDO
2023       ENDIF
2024 #else
2026         DO J=JDS,JDE-1
2027           IENDX=IDE-1
2028           IF(MOD(J,2)==0)IENDX=IENDX-1
2029           WRITE(IUNIT)(ARRAY(I,K,J),I=1,IENDX)
2030         ENDDO
2032 #endif
2034 !----------------------------------------------------------------------
2035   500 CONTINUE
2037       IF(MYPE==0)CLOSE(IUNIT)
2038 !----------------------------------------------------------------------
2040       END SUBROUTINE VWR
2041 !----------------------------------------------------------------------
2042 !**********************************************************************
2043 !----------------------------------------------------------------------
2044       SUBROUTINE EXIT(NAME,PINT,T,Q,U,V,Q2,W                           &
2045      &               ,NTSD,MYPE,MPI_COMM_COMP                          &
2046      &               ,IDS,IDE,JDS,JDE,KDS,KDE                          &
2047      &               ,IMS,IME,JMS,JME,KMS,KME                          &
2048      &               ,ITS,ITE,JTS,JTE,KTS,KTE)
2049 !----------------------------------------------------------------------
2050 !**********************************************************************
2051       USE MODULE_EXT_INTERNAL
2053 !----------------------------------------------------------------------
2054       IMPLICIT NONE
2055 !----------------------------------------------------------------------
2056 #if defined(DM_PARALLEL) && !defined(STUBMPI)
2057       INCLUDE "mpif.h"
2058 #endif
2059 !----------------------------------------------------------------------
2060       INTEGER,INTENT(IN) :: IDS,IDE,JDS,JDE,KDS,KDE                    &
2061      &                     ,IMS,IME,JMS,JME,KMS,KME                    &
2062      &                     ,ITS,ITE,JTS,JTE,KTS,KTE                    &
2063      &                     ,MYPE,MPI_COMM_COMP,NTSD
2065       REAL,DIMENSION(IMS:IME,JMS:JME,KMS:KME),INTENT(IN) :: PINT,T,Q   &
2066                                                            ,U,V,Q2,W
2067       CHARACTER(*),INTENT(IN) :: NAME
2069       INTEGER :: I,J,K,IEND,IERR,IRET
2070       CHARACTER(256) :: ERRMESS
2071       LOGICAL :: E_BDY,S_BDY
2072 !----------------------------------------------------------------------
2073       IRET=0
2074   100 FORMAT(' EXIT ',A,' AT NTSD=',I5)
2075       IEND=ITE
2076       S_BDY=(JTS==JDS)
2077       E_BDY=(ITE==IDE-1)
2079       DO K=KTS,KTE
2080       DO J=JTS,JTE
2081       IEND=ITE
2082       IF(E_BDY.AND.MOD(J,2)==0)IEND=ITE-1
2084       DO I=ITS,IEND
2085         IF(T(I,J,K)>330..OR.T(I,J,K)<180..OR.T(I,J,K)/=T(I,J,K))THEN
2086           WRITE(0,100)NAME,NTSD
2087           WRITE(0,200)I,J,K,T(I,J,K),MYPE,NTSD
2088   200     FORMAT(' BAD VALUE I=',I3,' J=',I3,' K=',I2,' T=',E12.5      &
2089      &,          ' MYPE=',I3,' NTSD=',I5)
2090           IRET=666
2091           return
2092 !         WRITE(ERRMESS,205)NAME,T(I,J,K),I,J,K,MYPE
2093   205     FORMAT(' EXIT ',A,' TEMPERATURE=',E12.5                      &
2094      &,          ' AT (',I3,',',I2,',',I3,')',' MYPE=',I3)
2095 !         CALL WRF_ERROR_FATAL(ERRMESS)
2096 !         CALL MPI_ABORT(MPI_COMM_WORLD,1,IERR)
2097         ELSEIF(Q(I,J,K)<-1.E-4.OR.Q(I,J,K)>30.E-3                      &
2098      &         .OR.Q(I,J,K)/=Q(I,J,K))THEN
2099           WRITE(0,100)NAME,NTSD
2100           WRITE(0,300)I,J,K,Q(I,J,K),MYPE,NTSD
2101   300     FORMAT(' BAD VALUE I=',I3,' J=',I3,' K=',I2,' Q=',E12.5      &
2102      &,          ' MYPE=',I3,' NTSD=',I5)
2103           IRET=666
2104           return
2105 !         WRITE(ERRMESS,305)NAME,Q(I,J,K),I,J,K,MYPE
2106   305     FORMAT(' EXIT ',A,' SPEC HUMIDITY=',E12.5                    &
2107      &,          ' AT (',I3,',',I2,',',I3,')',' MYPE=',I3)
2108 !         CALL WRF_ERROR_FATAL(ERRMESS)
2109 !         CALL MPI_ABORT(MPI_COMM_WORLD,1,IERR)
2110         ELSEIF(PINT(I,J,K)<0..OR.PINT(I,J,K)/=PINT(I,J,K))THEN
2111           WRITE(0,100)NAME,NTSD
2112           WRITE(0,315)I,J,K,PINT(I,J,K),MYPE,NTSD
2113   315     FORMAT(' BAD VALUE I=',I3,' J=',I3,' K=',I2,' PINT=',E12.5      &
2114      &,          ' MYPE=',I3,' NTSD=',I5)
2115           IRET=666
2116           return
2117 !         CALL WRF_ERROR_FATAL(ERRMESS)
2118 !         CALL MPI_ABORT(MPI_COMM_WORLD,1,IERR)
2119         ELSEIF(W(I,J,K)/=W(I,J,K))THEN
2120           WRITE(0,100)NAME,NTSD
2121           WRITE(0,325)I,J,K,W(I,J,K),MYPE,NTSD
2122   325     FORMAT(' BAD VALUE I=',I3,' J=',I3,' K=',I2,' W=',E12.5      &
2123      &,          ' MYPE=',I3,' NTSD=',I5)
2124           IRET=666
2125           return
2126 !         CALL WRF_ERROR_FATAL(ERRMESS)
2127 !         CALL MPI_ABORT(MPI_COMM_WORLD,1,IERR)
2128         ENDIF
2129       ENDDO
2130       ENDDO
2131       ENDDO
2133       DO K=KTS,KTE
2134       DO J=JTS,JTE
2135       IEND=ITE
2136       IF(E_BDY.AND.MOD(J,2)==1)IEND=ITE-1
2137       DO I=ITS,IEND
2138         IF(ABS(U(I,J,K))>125..OR.ABS(V(I,J,K))>125.                    &
2139      &         .OR.U(I,J,K)/=U(I,J,K).OR.V(I,J,K)/=V(I,J,K))THEN
2140           WRITE(0,100)NAME,NTSD
2141           WRITE(0,400)I,J,K,U(I,J,K),V(I,J,K),MYPE,NTSD
2142   400     FORMAT(' BAD VALUE I=',I3,' J=',I3,' K=',I2,' U=',E12.5      &
2143      &,          ' V=',E12.5,' MYPE=',I3,' NTSD=',I5)
2144           IRET=666
2145           return
2146 !         WRITE(ERRMESS,405)NAME,U(I,J,K),V(I,J,K),I,J,K,MYPE
2147   405     FORMAT(' EXIT ',A,' U=',E12.5,' V=',E12.5                    &
2148      &,          ' AT (',I3,',',I2,',',I3,')',' MYPE=',I3)
2149 !         CALL WRF_ERROR_FATAL(ERRMESS)
2150 !         CALL MPI_ABORT(MPI_COMM_WORLD,1,IERR)
2151         ENDIF
2152       ENDDO
2153       ENDDO
2154       ENDDO
2155 !----------------------------------------------------------------------
2156       END SUBROUTINE EXIT
2157 !----------------------------------------------------------------------
2158 !**********************************************************************
2159 !----------------------------------------------------------------------
2160       SUBROUTINE TIME_STATS(TIME_LCL,NAME,NTSD,MYPE,NPES,MPI_COMM_COMP)
2161 !----------------------------------------------------------------------
2162 !**********************************************************************
2163       USE MODULE_EXT_INTERNAL
2165 !----------------------------------------------------------------------
2166       IMPLICIT NONE
2167 !----------------------------------------------------------------------
2168 #if defined(DM_PARALLEL) && !defined(STUBMPI)
2169       INCLUDE "mpif.h"
2170 #endif
2171 !----------------------------------------------------------------------
2172       INTEGER,INTENT(IN) :: MPI_COMM_COMP,MYPE,NPES,NTSD
2173       REAL,INTENT(IN) :: TIME_LCL
2175       CHARACTER(*),INTENT(IN) :: NAME
2177 !*** LOCAL VARIABLES
2179 #if defined(DM_PARALLEL) && !defined(STUBMPI)
2180       INTEGER,DIMENSION(MPI_STATUS_SIZE) :: JSTAT
2181       INTEGER,DIMENSION(MPI_STATUS_SIZE,4) :: STATUS_ARRAY
2182 #endif
2183       INTEGER,ALLOCATABLE,DIMENSION(:) :: ID_PE,IPE_SORT
2185       INTEGER :: IPE,IPE_MAX,IPE_MEDIAN,IPE_MIN,IRECV,IRTN,ISEND       &
2186      &          ,N,N_MEDIAN,NLEN
2188       REAL,ALLOCATABLE,DIMENSION(:) :: TIME,SORT_TIME
2189       REAL,DIMENSION(2) :: REMOTE
2190       REAL :: TIME_MAX,TIME_MEAN,TIME_MEDIAN,TIME_MIN
2192       CHARACTER(5) :: TIMESTEP
2193       CHARACTER(6) :: FMT
2194       CHARACTER(25) :: TITLE
2195 !----------------------------------------------------------------------
2196 !**********************************************************************
2197 !----------------------------------------------------------------------
2199       IF(NTSD<=9)THEN
2200         FMT='(I1.1)'
2201         NLEN=1
2202       ELSEIF(NTSD<=99)THEN
2203         FMT='(I2.2)'
2204         NLEN=2
2205       ELSEIF(NTSD<=999)THEN
2206         FMT='(I3.3)'
2207         NLEN=3
2208       ELSEIF(NTSD<=9999)THEN
2209         FMT='(I4.4)'
2210         NLEN=4
2211       ELSEIF(NTSD<=99999)THEN
2212         FMT='(I5.5)'
2213         NLEN=5
2214       ENDIF
2215       WRITE(TIMESTEP,FMT)NTSD
2216       TITLE=NAME//'_'//TIMESTEP(1:NLEN)
2218 !----------------------------------------------------------------------
2220 #if defined(DM_PARALLEL) && !defined(STUBMPI)
2221       IF(MYPE==0)THEN
2222         ALLOCATE(TIME(1:NPES))
2223         ALLOCATE(SORT_TIME(1:NPES))
2224         ALLOCATE(ID_PE(1:NPES))
2225         ALLOCATE(IPE_SORT(1:NPES))
2227         TIME(1)=TIME_LCL
2228         ID_PE(1)=MYPE
2230 !***  COLLECT TIMES AND PE VALUES FROM OTHER PEs
2232         DO IPE=1,NPES-1
2233           CALL MPI_RECV(REMOTE,2,MPI_REAL,IPE,IPE                      &
2234      &                 ,MPI_COMM_COMP,JSTAT,IRECV)
2236           TIME(IPE+1)=REMOTE(1)
2237           ID_PE(IPE+1)=NINT(REMOTE(2))
2238         ENDDO
2240 !***  NOW GET STATS.
2241 !***  FIRST THE MAX, MIN, AND MEAN TIMES.
2243         TIME_MEAN=0.
2244         TIME_MAX=-1.
2245         TIME_MIN=1.E10
2246         IPE_MAX=-1
2247         IPE_MIN=-1
2249         DO N=1,NPES
2250           TIME_MEAN=TIME_MEAN+TIME(N)
2252           IF(TIME(N)>TIME_MAX)THEN
2253             TIME_MAX=TIME(N)
2254             IPE_MAX=ID_PE(N)
2255           ENDIF
2257           IF(TIME(N)<TIME_MIN)THEN
2258             TIME_MIN=TIME(N)
2259             IPE_MIN=ID_PE(N)
2260           ENDIF
2262         ENDDO
2264         TIME_MAX=TIME_MAX*1.E-3
2265         TIME_MIN=TIME_MIN*1.E-3
2266         TIME_MEAN=TIME_MEAN*1.E-3/REAL(NPES)
2268 !***  THEN THE MEDIAN TIME.
2270         CALL SORT(TIME,NPES,SORT_TIME,IPE_SORT)
2271         N_MEDIAN=(NPES+1)/2
2272         TIME_MEDIAN=SORT_TIME(N_MEDIAN)*1.E-3
2273         IPE_MEDIAN=IPE_SORT(N_MEDIAN)
2275 !----------------------------------------------------------------------
2276       ELSE
2278 !***  SEND TIME AND PE VALUE TO PE0.
2280         REMOTE(1)=TIME_LCL
2281         REMOTE(2)=REAL(MYPE)
2283         CALL MPI_SEND(REMOTE,2,MPI_REAL,0,MYPE,MPI_COMM_COMP,ISEND)
2285       ENDIF
2286 !----------------------------------------------------------------------
2288       CALL MPI_BARRIER(MPI_COMM_COMP,IRTN)
2290 !***  WRITE RESULTS
2292       IF(MYPE==0)THEN
2293         WRITE(0,100)TITLE
2294         WRITE(0,105)TIME_MAX,IPE_MAX
2295         WRITE(0,110)TIME_MIN,IPE_MIN
2296         WRITE(0,115)TIME_MEDIAN,IPE_MEDIAN
2297         WRITE(0,120)TIME_MEAN
2298   100   FORMAT(' Time for ',A)
2299   105   FORMAT(' Maximum=',G11.5,' for PE ',I2.2)
2300   110   FORMAT(' Minimum=',G11.5,' for PE ',I2.2)
2301   115   FORMAT(' Median =',G11.5,' for PE ',I2.2)
2302   120   FORMAT(' Mean   =',G11.5)
2303       ENDIF
2304 !----------------------------------------------------------------------
2306 #endif
2307       END SUBROUTINE TIME_STATS
2309 !----------------------------------------------------------------------
2310 !**********************************************************************
2311 !----------------------------------------------------------------------
2312       SUBROUTINE SORT(DATA,NPES,DATA_SORTED,IPE_SORTED)
2313 !----------------------------------------------------------------------
2314 !***
2315 !***  SORT DATA FROM MULTIPLE PEs.  SEND BACK THE SORTED DATA ITEMS
2316 !***  ALONG WITH THE ASSOCIATED TASK IDs.
2317 !***
2318 !----------------------------------------------------------------------
2319       IMPLICIT NONE
2320 !----------------------------------------------------------------------
2321       INTEGER,INTENT(IN) :: NPES
2322       REAL,DIMENSION(NPES),INTENT(IN) :: DATA
2324       INTEGER,DIMENSION(NPES),INTENT(OUT) :: IPE_SORTED
2325       REAL,DIMENSION(NPES),INTENT(OUT) :: DATA_SORTED
2326 !----------------------------------------------------------------------
2327       TYPE :: DATA_LINK
2328         REAL :: VALUE
2329         INTEGER :: IPE
2330         TYPE(DATA_LINK),POINTER :: NEXT_VALUE
2331       END TYPE
2332 !----------------------------------------------------------------------
2334 !***  LOCAL VARIABLES
2336 !----------------------------------------------------------------------
2337       INTEGER :: ISTAT,N
2339       TYPE(DATA_LINK),POINTER :: HEAD,TAIL  ! Smallest, largest
2340       TYPE(DATA_LINK),POINTER :: PTR_NEW    ! Each new value
2341       TYPE(DATA_LINK),POINTER :: PTR1,PTR2  ! Working pointers
2342 !----------------------------------------------------------------------
2343 !**********************************************************************
2344 !----------------------------------------------------------------------
2345       pe_loop: DO N=1,NPES
2346         ALLOCATE(PTR_NEW,STAT=ISTAT)  ! Location for next data items
2347         PTR_NEW%VALUE=DATA(N)
2348         PTR_NEW%IPE=N-1
2350 !----------------------------------------------------------------------
2351 !***  DETERMINE WHERE IN LIST TO INSERT VALUE.
2352 !***  FIRST THE INITIAL DATA VALUE.
2353 !----------------------------------------------------------------------
2355 !       main: IF(.NOT.ASSOCIATED(HEAD))THEN
2356         main: IF(N==1)THEN
2357           HEAD=>PTR_NEW
2358           TAIL=>HEAD
2359           NULLIFY(PTR_NEW%NEXT_VALUE)
2361 !----------------------------------------------------------------------
2362 !***  THE NEW VALUE IS LESS THAN THE SMALLEST VALUE ALREADY SORTED.
2363 !----------------------------------------------------------------------
2365         ELSE
2366           check: IF(PTR_NEW%VALUE<HEAD%VALUE)THEN
2367             PTR_NEW%NEXT_VALUE=>HEAD
2368             HEAD=>PTR_NEW
2370 !----------------------------------------------------------------------
2371 !***  THE NEW VALUE IS GREATER THAN THE LARGEST VALUE ALREADY SORTED.
2372 !----------------------------------------------------------------------
2374           ELSEIF(PTR_NEW%VALUE>=TAIL%VALUE)THEN
2375             TAIL%NEXT_VALUE=>PTR_NEW  ! This is what connects the former
2376                                       ! final value in the list to
2377                                       ! the new value being appended.
2378             TAIL=>PTR_NEW
2379             NULLIFY(TAIL%NEXT_VALUE)
2381 !----------------------------------------------------------------------
2382 !***  THE NEW VALUE IS IN BETWEEN VALUES ALREADY SORTED.
2383 !----------------------------------------------------------------------
2385           ELSE
2386             PTR1=>HEAD
2387             PTR2=>PTR1%NEXT_VALUE
2389             search: DO
2390               IF((PTR_NEW%VALUE>=PTR1%VALUE).AND.                      &
2391      &           (PTR_NEW%VALUE<PTR2%VALUE))THEN
2392                 PTR_NEW%NEXT_VALUE=>PTR2
2393                 PTR1%NEXT_VALUE=>PTR_NEW
2394                 EXIT search
2395               ENDIF
2397               PTR1=>PTR2
2398               PTR2=>PTR2%NEXT_VALUE
2399             ENDDO search
2401           ENDIF check
2403         ENDIF main
2405       ENDDO pe_loop
2407 !----------------------------------------------------------------------
2408 !***  COLLECT THE SORTED NUMBERS FROM THE LINKED LIST.
2409 !----------------------------------------------------------------------
2411       PTR1=>HEAD
2413       DO N=1,NPES
2414 !       IF(.NOT.ASSOCIATED(PTR_NEW))EXIT
2415         DATA_SORTED(N)=PTR1%VALUE
2416         IPE_SORTED(N)=PTR1%IPE
2417         PTR1=>PTR1%NEXT_VALUE
2418       ENDDO
2420       DEALLOCATE(PTR_NEW)
2421       NULLIFY (HEAD,TAIL,PTR1,PTR2)
2422 !----------------------------------------------------------------------
2423       END SUBROUTINE SORT
2424 !----------------------------------------------------------------------
2425 !**********************************************************************
2426 !----------------------------------------------------------------------
2427       SUBROUTINE FIELD_STATS(FIELD,MYPE,MPI_COMM_COMP                  &
2428      &                      ,IDS,IDE,JDS,JDE,KDS,KDE                   &
2429      &                      ,IMS,IME,JMS,JME,KMS,KME                   &
2430      &                      ,ITS,ITE,JTS,JTE,KTS,KTE)
2431 !----------------------------------------------------------------------
2432 !***
2433 !***  GENERATE STANDARD LAYER STATISTICS FOR THE DESIRED FIELD.
2434 !***
2435 !----------------------------------------------------------------------
2436       IMPLICIT NONE
2437 !----------------------------------------------------------------------
2438 #if defined(DM_PARALLEL) && !defined(STUBMPI)
2439       INCLUDE "mpif.h"
2440 #endif
2441 !----------------------------------------------------------------------
2443       INTEGER,INTENT(IN) :: MPI_COMM_COMP,MYPE
2444       INTEGER,INTENT(IN) :: IDS,IDE,JDS,JDE,KDS,KDE                    &
2445      &                     ,IMS,IME,JMS,JME,KMS,KME                    &
2446      &                     ,ITS,ITE,JTS,JTE,KTS,KTE
2448       REAL,DIMENSION(IMS:IME,JMS:JME,KMS:KME),INTENT(IN) :: FIELD
2450 !----------------------------------------------------------------------
2451 !***  LOCAL
2452 !----------------------------------------------------------------------
2454       INTEGER,PARAMETER :: DOUBLE=SELECTED_REAL_KIND(15,300)
2456       INTEGER :: I,IEND,IRTN,I_BY_J,J,K,KFLIP
2458       REAL :: FIKJ,FMAXK,FMINK
2459       REAL(KIND=DOUBLE) :: F_MEAN,POINTS,RMS,ST_DEV,SUMFK,SUMF2K
2460       REAL,DIMENSION(KTS:KTE) :: FMAX,FMAX_0,FMIN,FMIN_0
2461       REAL(KIND=DOUBLE),DIMENSION(KTS:KTE) :: SUMF,SUMF_0,SUMF2,SUMF2_0
2462 !----------------------------------------------------------------------
2464       I_BY_J=(IDE-IDS)*(JDE-JDS)-(JDE-JDS-1)/2  ! This assumes that
2465                                                 ! IDE and JDE are each
2466                                                 ! one greater than the
2467                                                 ! true grid size.
2468 #if defined(DM_PARALLEL) && !defined(STUBMPI)
2470       layer_loop:  DO K=KTS,KTE
2472         FMAXK=-1.E10
2473         FMINK=1.E10
2474         SUMFK=0.
2475         SUMF2K=0.
2477         DO J=JTS,JTE
2478           IEND=MIN(ITE,IDE-1)
2479           IF(MOD(J,2)==0.AND.ITE==IDE-1)IEND=IEND-1
2480           DO I=ITS,IEND
2481             FIKJ=FIELD(I,J,K)
2482             FMAXK=MAX(FMAXK,FIKJ)
2483             FMINK=MIN(FMINK,FIKJ)
2484             SUMFK=SUMFK+FIKJ
2485             SUMF2K=SUMF2K+FIKJ*FIKJ
2486           ENDDO
2487         ENDDO
2489         FMAX(K)=FMAXK
2490         FMIN(K)=FMINK
2491         SUMF(K)=SUMFK
2492         SUMF2(K)=SUMF2K
2494       ENDDO layer_loop
2496 !----------------------------------------------------------------------
2497 !***  GLOBAL STATS
2498 !----------------------------------------------------------------------
2500       CALL MPI_REDUCE(SUMF,SUMF_0,KTE,MPI_REAL8,MPI_SUM,0              &
2501      &               ,MPI_COMM_COMP,IRTN)
2502       CALL MPI_REDUCE(SUMF2,SUMF2_0,KTE,MPI_REAL8,MPI_SUM,0            &
2503      &               ,MPI_COMM_COMP,IRTN)
2504       CALL MPI_REDUCE(FMAX,FMAX_0,KTE,MPI_REAL,MPI_MAX,0               &
2505      &               ,MPI_COMM_COMP,IRTN)
2506       CALL MPI_REDUCE(FMIN,FMIN_0,KTE,MPI_REAL,MPI_MIN,0               &
2507      &               ,MPI_COMM_COMP,IRTN)
2509       IF(MYPE==0)THEN
2510         POINTS=I_BY_J
2511         DO K=KTE,KTS,-1
2512           F_MEAN=SUMF_0(K)/POINTS
2513           ST_DEV=SQRT((POINTS*SUMF2_0(K)-SUMF_0(K)*SUMF_0(K))/         &
2514      &                (POINTS*(POINTS-1)))
2515           RMS=SQRT(SUMF2_0(K)/POINTS)
2516           KFLIP=KTE-K+1
2517           WRITE(0,101)KFLIP,FMAX_0(K),FMIN_0(K)
2518           WRITE(0,102)F_MEAN,ST_DEV,RMS
2519   101     FORMAT(' LAYER=',I2,' MAX=',E13.6,' MIN=',E13.6)
2520   102     FORMAT(9X,' MEAN=',E13.6,' STDEV=',E13.6,' RMS=',E13.6)
2521         ENDDO
2522       ENDIF
2523 #endif
2524 !----------------------------------------------------------------------
2525       END SUBROUTINE FIELD_STATS
2526 !----------------------------------------------------------------------
2527       FUNCTION TIMEF()
2528       REAL*8 TIMEF
2529       INTEGER :: IC,IR
2530       CALL SYSTEM_CLOCK(COUNT=IC,COUNT_RATE=IR)
2531       TIMEF=REAL(IC)/REAL(IR)*1000.0
2532       END