UPP  001
 All Data Structures Files Functions Pages
WRFPOST.f
Go to the documentation of this file.
1 
32  PROGRAM wrfpost
33 
34 !
35 !
36 !============================================================================================================
37 !
38 ! This is an MPI code. All array indexing is with respect to the global indices. Loop indices
39 ! look as follows for N MPI tasks.
40 !
41 !
42 !
43 ! Original New
44 ! Index Index
45 !
46 ! JM ----------------------------------------------- JEND
47 ! JM-1 - - JEND_M
48 ! JM-2 - MPI TASK N-1 - JEND_M2
49 ! - -
50 ! - -
51 ! ----------------------------------------------- JSTA, JSTA_M, JSTA_M2
52 ! ----------------------------------------------- JEND, JEND_M, JEND_M2
53 ! - -
54 ! - MPI TASK N-2 -
55 ! - -
56 ! - -
57 ! ----------------------------------------------- JSTA, JSTA_M, JSTA_M2
58 !
59 ! .
60 ! .
61 ! .
62 !
63 ! ----------------------------------------------- JEND, JEND_M, JEND_M2
64 ! - -
65 ! - MPI TASK 1 -
66 ! - -
67 ! - -
68 ! ----------------------------------------------- JSTA, JSTA_M, JSTA_M2
69 ! ----------------------------------------------- JEND, JEND_M, JEND_M2
70 ! - -
71 ! - MPI TASK 0 -
72 ! 3 - - JSTA_M2
73 ! 2 - - JSTA_M
74 ! 1 ----------------------------------------------- JSTA
75 !
76 ! 1 IM
77 !
78 !
79 ! Jim Tuccillo
80 ! Jan 2000
81 !
82 ! README - Jim Tuccillo Feb 2001
83 !
84 ! Many common blocks have been replaced by modules to support Fortran
85 ! "allocate" commands. Many of the 3-D arrays are now allocated to be the
86 ! exact size required based on the number of MPI tasks. The dimensioning will be
87 ! x ( im,jsta_2l:jend_2u,lm)
88 ! Most 2-D arrays continue to be dimensioned (im,jm). This is fine but please be aware
89 ! that the EXCH routine for arrays dimensioned (im,jm) is different than arrays dimensioned
90 ! (im,jsta_2l:jend_2u). Also, be careful about passing any arrays dimensioned
91 ! (im,jst_2l:jend_2u,lm). See examples in the code as to the correct calling sequence and
92 ! EXCH routine to use.
93 !
94 !
95 ! ASYNCHRONOUS I/O HAS BEEN ADDED. THE LAST MPI TASK DOES THE I/O. IF THERE IS
96 ! ONLY ONE MPI TASK THN TASK ) DOES THE I/O.
97 ! THE CODE HAS GOTTEN A LITTLE KLUDGY. BASICLY, IM, IMX and IMOUT MUST BE EQUAL
98 ! AND REPRESENT THE VALUE USED IN THE MODEL. THE SAME HOLDS FOR JM, JMX and JMOUT.
99 !
100 ! Jim Tuccillo June 2001
101 !
102 !
103 !===========================================================================================
104 !
105  use netcdf
106  use nemsio_module, only: nemsio_getheadvar, nemsio_gfile, nemsio_init, nemsio_open, &
107  nemsio_getfilehead,nemsio_close
108  use ctlblk_mod, only: filenameaer, me, num_procs, num_servers, mpi_comm_comp, datestr, &
109  mpi_comm_inter, filename, ioform, grib, idat, filenameflux, filenamed3d, gdsdegr, &
110  spldef, modelname, ihrst, lsmdef,vtimeunits, tprec, pthresh, datahandle, im, jm, lm, &
111  lp1, lm1, im_jm, isf_surface_physics, nsoil, spl, lsmp1, global, imp_physics, &
112  jsta, jend, jsta_m, jend_m, jsta_2l, jend_2u, novegtype, icount_calmict, npset, datapd,&
113  lsm, fld_info, etafld2_tim, eta2p_tim, mdl2sigma_tim, cldrad_tim, miscln_tim, &
114  mdl2agl_tim, mdl2std_tim, mdl2thandpv_tim, calrad_wcloud_tim, &
115  fixed_tim, time_output, imin, surfce2_tim, komax, ivegsrc, d3d_on, gocart_on,rdaod, &
116  readxml_tim, spval, fullmodelname, submodelname, hyb_sigp, filenameflat, aqfcmaq_on
117  use grib2_module, only: gribit2,num_pset,nrecout,first_grbtbl,grib_info_finalize
118 !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
119  implicit none
120 !
121  type(nemsio_gfile) :: nfile,ffile,rfile
122  include "mpif.h"
123 !
124 ! DECLARE VARIABLES.
125 !
126 ! SET HEADER WRITER FLAGS TO TRUE.
127 !
128 !temporary vars
129 !
130  real(kind=8) :: time_initpost=0.,initpost_tim=0.,btim,bbtim
131  real rinc(5), untcnvt
132  integer :: status=0,iostatusd3d=0,iostatusflux=0
133  integer i,j,iii,l,k,ierr,nrec,ist,lusig,idrt,ncid3d,ncid2d,varid
134  integer :: prntsec,iim,jjm,llm,ioutcount,itmp,iret,iunit, &
135  iunitd3d,iyear,imn,iday,lcntrl,ieof
136  integer :: iostatusaer
137  logical :: popascal
138 !
139  integer :: kpo,kth,kpv
140  real,dimension(komax) :: po,th,pv
141  namelist/nampgb/kpo,po,kth,th,kpv,pv,filenameaer,d3d_on,gocart_on,popascal &
142  ,hyb_sigp,rdaod,aqfcmaq_on,vtimeunits
143  integer :: itag_ierr
144  namelist/model_inputs/filename,ioform,grib,datestr,modelname,submodelname &
145  ,filenameflux,filenameflat
146 
147  character startdate*19,sysdepinfo*80,iowrfname*3,post_fname*255
148  character cgar*1,cdum*4,line*10
149 !
150 !------------------------------------------------------------------------------
151 ! START HERE
152 !
153  call start()
154 !
155 ! INITIALIZE MPI
156 
157  CALL setup_servers(me, &
158  & num_procs, &
159  & num_servers, &
160  & mpi_comm_comp, &
161  & mpi_comm_inter)
162 !
163 ! ME IS THE RANK
164 ! NUM_PROCS IS THE NUMBER OF TASKS DOING POSTING
165 ! NUM_SERVERS IS ONE IF THERE ARE MORE THAN ONE TOTAL MPI TASKS, OTHERWISE ZERO
166 ! MPI_COMM_COMP IS THE INTRACOMMUNICATOR
167 ! MPI_COMM_INTER IS THE INTERCOMMUNICATOR FOR COMMUNCATION BETWEEN TASK 0 OF THE
168 ! TASKS DOING THE POSTING AND THE I/O SERVER
169 !
170 !
171 ! IF WE HAVE MORE THAN 1 MPI TASK THEN WE WILL FIRE UP THE IO SERVER
172 ! THE LAST TASK ( IN THE CONTEXT OF MPI_COMM_WORLD ) IS THE I/O SERVER
173 !
174  print*,'ME,NUM_PROCS,NUM_SERVERS=',me,num_procs,num_servers
175 
176  if (me == 0) CALL w3tagb('nems ',0000,0000,0000,'np23 ')
177 
178  if ( me >= num_procs ) then
179 !
180  call server
181 !
182  else
183  spval = 9.9e10
184 !
185 !**************************************************************************
186 !KaYee: Read itag in Fortran Namelist format
187 !Set default
188  submodelname='NONE'
189 !open namelist
190  open(5,file='itag')
191  read(5,nml=model_inputs,iostat=itag_ierr,err=888)
192  !print*,'itag_ierr=',itag_ierr
193 888 if (itag_ierr /= 0) then
194  print*,'Incorrect namelist variable(s) found in the itag file,stopping!'
195  stop
196  endif
197  if (me==0) print*,'fileName= ',filename
198  if (me==0) print*,'IOFORM= ',ioform
199  !if (me==0) print*,'OUTFORM= ',grib
200  if (me==0) print*,'OUTFORM= ',grib
201  if (me==0) print*,'DateStr= ',datestr
202  if (me==0) print*,'MODELNAME= ',modelname
203  if (me==0) print*,'SUBMODELNAME= ',submodelname
204 ! if(MODELNAME == 'NMM')then
205 ! read(5,1114) VTIMEUNITS
206 ! 1114 format(a4)
207 ! if (me==0) print*,'VALID TIME UNITS = ', VTIMEUNITS
208 ! endif
209 !
210  303 format('MODELNAME="',a,'" SUBMODELNAME="',a,'"')
211 
212  write(0,*)'MODELNAME: ', modelname, submodelname
213 
214  if (me==0) print 303,modelname,submodelname
215 ! assume for now that the first date in the stdin file is the start date
216  read(datestr,300) iyear,imn,iday,ihrst,imin
217  if (me==0) write(*,*) 'in WRFPOST iyear,imn,iday,ihrst,imin', &
218  iyear,imn,iday,ihrst,imin
219  300 format(i4,1x,i2,1x,i2,1x,i2,1x,i2)
220 
221  idat(1) = imn
222  idat(2) = iday
223  idat(3) = iyear
224  idat(4) = ihrst
225  idat(5) = imin
226 
227  111 format(a256)
228  112 format(a19)
229  113 format(a20)
230  114 format(a8)
231  120 format(a5)
232  121 format(a4)
233 
234 !KaYee: Read in GFS/FV3 runs in Fortran Namelist Format.
235  if (me==0) print*,'MODELNAME= ',modelname,'grib=',grib
236  if(modelname == 'GFS' .OR. modelname == 'FV3R') then
237  if (me == 0) print*,'first two file names in GFS or FV3= ' &
238  ,trim(filename),trim(filenameflux)
239  end if
240 
241  if(grib=='grib2') then
242  gdsdegr = 1.d6
243  endif
244  if (me==0) print *,'gdsdegr=',gdsdegr
245 !
246 ! set default for kpo, kth, th, kpv, pv
247  kpo = 0
248  po = 0
249  kth = 6
250  th = (/310.,320.,350.,450.,550.,650.,(0.,k=kth+1,komax)/) ! isentropic level to output
251  kpv = 8
252  pv = (/0.5,-0.5,1.0,-1.0,1.5,-1.5,2.0,-2.0,(0.,k=kpv+1,komax)/)
253 
254  hyb_sigp = .true.
255  d3d_on = .false.
256  gocart_on = .false.
257  aqfcmaq_on = .false.
258  popascal = .false.
259  filenameaer = ''
260  rdaod = .false.
261 ! gocart_on = .true.
262 ! d3d_on = .true.
263 
264 !set control file name
265  filenameflat='postxconfig-NT.txt'
266  read(5,nampgb,iostat=iret,end=119)
267  119 continue
268  if(me == 0) then
269  print*,'komax,iret for nampgb= ',komax,iret
270  print*,'komax,kpo,kth,th,kpv,pv,fileNameAER,popascal= ',komax,kpo &
271  & ,kth,th(1:kth),kpv,pv(1:kpv),trim(filenameaer),popascal
272  endif
273 
274 ! set up pressure level from POSTGPVARS or DEFAULT
275  if(kpo == 0) then
276 ! use default pressure levels
277  if(me == 0) then
278  print*,'using default pressure levels,spldef=',(spldef(l),l=1,lsmdef)
279  endif
280  lsm = lsmdef
281  do l=1,lsm
282  spl(l) = spldef(l)
283  end do
284  else
285 ! use POSTGPVARS
286  if(me == 0) then
287  print*,'using pressure levels from POSTGPVARS'
288  endif
289  lsm = kpo
290  if( .not. popascal ) then
291  untcnvt = 100.
292  else
293  untcnvt = 1.
294  endif
295  if(po(lsm) < po(1))then ! post logic assumes asscending
296  do l=1,lsm
297  spl(l) = po(lsm-l+1)*untcnvt
298  end do
299  else
300  do l=1,lsm
301  spl(l) = po(l)*untcnvt
302  end do
303  end if
304  end if
305  lsmp1 = lsm+1
306  if (me==0) print*,'LSM, SPL = ',lsm,spl(1:lsm)
307 
308  116 continue
309 
310 ! set PTHRESH for different models
311  if(modelname == 'NMM')then
312  pthresh = 0.000004
313  else
314  pthresh = 0.000001
315  end if
316 !Chuang: add dynamical allocation
317  if(trim(ioform) == 'netcdf' .OR. trim(ioform) == 'netcdfpara') THEN
318  IF(modelname == 'NCAR' .OR. modelname == 'RAPR' .OR. modelname == 'NMM') THEN
319  call ext_ncd_ioinit(sysdepinfo,status)
320  print*,'called ioinit', status
321  call ext_ncd_open_for_read( trim(filename), 0, 0, " ", &
322  datahandle, status)
323  print*,'called open for read', status
324  if ( status /= 0 ) then
325  print*,'error opening ',filename, ' Status = ', status ; stop
326  endif
327  call ext_ncd_get_dom_ti_integer(datahandle &
328  ,'WEST-EAST_GRID_DIMENSION',iim,1,ioutcount, status )
329  im = iim-1
330  call ext_ncd_get_dom_ti_integer(datahandle &
331  ,'SOUTH-NORTH_GRID_DIMENSION',jjm,1,ioutcount, status )
332  jm = jjm-1
333  call ext_ncd_get_dom_ti_integer(datahandle &
334  ,'BOTTOM-TOP_GRID_DIMENSION',llm,1,ioutcount, status )
335  lm = llm-1
336  lp1 = lm+1
337  lm1 = lm-1
338  im_jm = im*jm
339 
340  print*,'im jm lm from wrfout= ',im,jm, lm
341 
342 ! Read and set global value for surface physics scheme
343  call ext_ncd_get_dom_ti_integer(datahandle &
344  ,'SF_SURFACE_PHYSICS',itmp,1,ioutcount, status )
345  isf_surface_physics = itmp
346  print*,'SF_SURFACE_PHYSICS= ',isf_surface_physics
347 ! set NSOIL to 4 as default for NOAH but change if using other
348 ! SFC scheme
349  nsoil = 4
350  IF(itmp == 1) then !thermal diffusion scheme
351  nsoil = 5
352  ELSE IF(itmp == 3) then ! RUC LSM
353  nsoil = 9
354  ELSE IF(itmp == 7) then ! Pleim Xu
355  nsoil = 2
356  END IF
357  print*,'NSOIL from wrfout= ',nsoil
358 
359  call ext_ncd_ioclose( datahandle, status )
360  ELSE
361 ! use parallel netcdf lib directly to read FV3 output in netCDF
362  spval = 9.99e20
363  status = nf90_open(trim(filename),ior(nf90_nowrite,nf90_mpiio), &
364  ncid3d,comm=mpi_comm_world,info=mpi_info_null)
365  if ( status /= 0 ) then
366  print*,'error opening ',filename, ' Status = ', status
367  stop
368  endif
369  status = nf90_open(trim(filenameflux),ior(nf90_nowrite,nf90_mpiio), &
370  ncid2d,comm=mpi_comm_world,info=mpi_info_null)
371  if ( status /= 0 ) then
372  print*,'error opening ',filenameflux, ' Status = ', status
373  stop
374  endif
375 ! read in LSM index and nsoil here
376  status=nf90_get_att(ncid2d,nf90_global,'landsfcmdl', isf_surface_physics)
377  if(status/=0)then
378  print*,'landsfcmdl not found; assigning to 2'
379  isf_surface_physics=2 !set LSM physics to 2 for NOAH
380  endif
381  if(isf_surface_physics<2)then
382  isf_surface_physics=2 !set LSM physics to 2 for NOAH
383  endif
384  status=nf90_get_att(ncid2d,nf90_global,'nsoil', nsoil)
385  if(status/=0)then
386  print*,'nsoil not found; assigning to 4'
387  nsoil=4 !set nsoil to 4 for NOAH
388  endif
389  if(me==0)print*,'SF_SURFACE_PHYSICS= ',isf_surface_physics
390  if(me==0)print*,'NSOIL= ',nsoil
391 ! read imp_physics
392  status=nf90_get_att(ncid2d,nf90_global,'imp_physics',imp_physics)
393  if(status/=0)then
394  print*,'imp_physics not found; assigning to GFDL 11'
395  imp_physics=11
396  endif
397  if (me == 0) print*,'MP_PHYSICS= ',imp_physics
398 ! get dimesions
399  status = nf90_inq_dimid(ncid3d,'grid_xt',varid)
400  if ( status /= 0 ) then
401  print*,status,varid
402  stop 1
403  end if
404  status = nf90_inquire_dimension(ncid3d,varid,len=im)
405  if ( status /= 0 ) then
406  print*,status
407  stop 1
408  end if
409  status = nf90_inq_dimid(ncid3d,'grid_yt',varid)
410  if ( status /= 0 ) then
411  print*,status,varid
412  stop 1
413  end if
414  status = nf90_inquire_dimension(ncid3d,varid,len=jm)
415  if ( status /= 0 ) then
416  print*,status
417  stop 1
418  end if
419  status = nf90_inq_dimid(ncid3d,'pfull',varid)
420  if ( status /= 0 ) then
421  print*,status,varid
422  stop 1
423  end if
424  status = nf90_inquire_dimension(ncid3d,varid,len=lm)
425  if ( status /= 0 ) then
426  print*,status
427  stop 1
428  end if
429  lp1 = lm+1
430  lm1 = lm-1
431  im_jm = im*jm
432 ! set NSOIL to 4 as default for NOAH but change if using other
433 ! SFC scheme
434 ! NSOIL = 4
435 
436  print*,'im jm lm nsoil from fv3 output = ',im,jm,lm,nsoil
437  END IF
438 
439  ELSE IF(trim(ioform) == 'binary' .OR. &
440  trim(ioform) == 'binarympiio' ) THEN
441  print*,'WRF Binary format is no longer supported'
442  stop 9996
443 ! NEMSIO format
444  ELSE IF(trim(ioform) == 'binarynemsio' .or. &
445  trim(ioform) == 'binarynemsiompiio' )THEN
446 
447  spval = 9.99e20
448  IF(me == 0)THEN
449  call nemsio_init(iret=status)
450  print *,'nemsio_init, iret=',status
451  call nemsio_open(nfile,trim(filename),'read',iret=status)
452  if ( status /= 0 ) then
453  print*,'error opening ',filename, ' Status = ', status ; stop
454  endif
455 !---
456  call nemsio_getfilehead(nfile,iret=status,nrec=nrec &
457  ,dimx=im,dimy=jm,dimz=lm,nsoil=nsoil)
458  if ( status /= 0 ) then
459  print*,'error finding model dimensions '; stop
460  endif
461  call nemsio_getheadvar(nfile,'global',global,iret)
462  if (iret /= 0)then
463  print*,"global not found in file-Assigned false"
464  global = .false.
465  end if
466  IF(modelname == 'GFS') global = .true.
467 ! global NMMB has i=1 overlap with i=im so post will cut off i=im
468  if(global .and. modelname == 'NMM') im = im-1
469 
470  END IF
471 
472  CALL mpi_bcast(im, 1,mpi_integer,0, mpi_comm_comp,status)
473  call mpi_bcast(jm, 1,mpi_integer,0, mpi_comm_comp,status)
474  call mpi_bcast(lm, 1,mpi_integer,0, mpi_comm_comp,status)
475  call mpi_bcast(nsoil,1,mpi_integer,0, mpi_comm_comp,status)
476 
477  if (me == 0) print*,'im jm lm nsoil from NEMS= ',im,jm, lm ,nsoil
478  call mpi_bcast(global,1,mpi_logical,0,mpi_comm_comp,status)
479  if (me == 0) print*,'Is this a global run ',global
480  lp1 = lm+1
481  lm1 = lm-1
482  im_jm = im*jm
483 
484 ! opening GFS flux file
485  IF(modelname == 'GFS') THEN
486 ! iunit=33
487  call nemsio_open(ffile,trim(filenameflux),'read',iret=iostatusflux)
488  if ( iostatusflux /= 0 ) then
489  print*,'error opening ',filenameflux, ' Status = ', iostatusflux
490  endif
491  iostatusd3d = -1
492  iunitd3d = -1
493 !
494 ! opening GFS aer file
495  call nemsio_open(rfile,trim(filenameaer),'read',iret=iostatusaer)
496  if ( iostatusaer /= 0 .and. me == 0) then
497  print*,'error opening AER ',filenameaer, ' Status = ', iostatusaer
498  endif
499 !
500 ! print*,'iostatusD3D in WRFPOST= ',iostatusD3D
501 
502  END IF
503 
504  ELSE
505  print*,'UNKNOWN MODEL OUTPUT FORMAT, STOPPING'
506  stop 9999
507  END IF
508 
509 
510  CALL mpi_first()
511  print*,'jsta,jend,jsta_m,jend_m,jsta_2l,jend_2u,spval=',jsta, &
512  jend,jsta_m,jend_m, jsta_2l,jend_2u,spval
513  CALL allocate_all()
514 
515 !
516 ! INITIALIZE POST COMMON BLOCKS
517 !
518  lcntrl = 14
519  rewind(lcntrl)
520 
521 ! EXP. initialize netcdf here instead
522  bbtim = mpi_wtime()
523  btim = mpi_wtime()
524 ! set default novegtype
525  if(modelname == 'GFS')THEN
526  novegtype = 13
527  ivegsrc = 2
528  else if(modelname=='NMM' .and. trim(ioform)=='binarynemsio')then
529  novegtype = 20
530  ivegsrc = 1
531  else if(modelname=='RAPR')then
532  novegtype = 20
533  ivegsrc = 1
534  else ! USGS
535  novegtype = 24
536  ivegsrc = 0
537  end if
538 
539 ! Reading model output for different models and IO format
540 
541  IF(trim(ioform) == 'netcdf' .OR. trim(ioform) == 'netcdfpara') THEN
542  IF(modelname == 'NCAR' .OR. modelname == 'RAPR') THEN
543  print*,'CALLING INITPOST TO PROCESS NCAR NETCDF OUTPUT'
544  CALL initpost
545  ELSE IF (modelname == 'FV3R' .OR. modelname == 'GFS') THEN
546 ! use parallel netcdf library to read output directly
547  print*,'CALLING INITPOST_NETCDF'
548  CALL initpost_netcdf(ncid2d,ncid3d)
549  ELSE
550  print*,'POST does not have netcdf option for model,',modelname,' STOPPING,'
551  stop 9998
552  END IF
553  ELSE IF(trim(ioform) == 'binarympiio') THEN
554  IF(modelname == 'NCAR' .OR. modelname == 'RAPR' .OR. modelname == 'NMM') THEN
555  print*,'WRF BINARY IO FORMAT IS NO LONGER SUPPORTED, STOPPING'
556  stop 9996
557  ELSE IF(modelname == 'RSM') THEN
558  print*,'MPI BINARY IO IS NOT YET INSTALLED FOR RSM, STOPPING'
559  stop 9997
560  ELSE
561  print*,'POST does not have mpiio option for this model, STOPPING'
562  stop 9998
563  END IF
564  ELSE IF(trim(ioform) == 'binarynemsio') THEN
565  IF(modelname == 'NMM') THEN
566  CALL initpost_nems(nrec,nfile)
567  ELSE
568  print*,'POST does not have nemsio option for model,',modelname,' STOPPING,'
569  stop 9998
570 
571  END IF
572 
573  ELSE IF(trim(ioform) == 'binarynemsiompiio')THEN
574  IF(modelname == 'GFS') THEN
575 ! close nemsio file for serial read
576  call nemsio_close(nfile,iret=status)
577  call nemsio_close(ffile,iret=status)
578  call nemsio_close(rfile,iret=status)
579  CALL initpost_gfs_nems_mpiio(iostatusaer)
580  ELSE
581  print*,'POST does not have nemsio mpi option for model,',modelname, &
582  'STOPPING,'
583  stop 9999
584 
585  END IF
586 
587  ELSE
588  print*,'UNKNOWN MODEL OUTPUT FORMAT, STOPPING'
589  stop 9999
590  END IF
591  initpost_tim = initpost_tim +(mpi_wtime() - btim)
592  IF(me == 0)THEN
593  WRITE(6,*)'WRFPOST: INITIALIZED POST COMMON BLOCKS'
594  ENDIF
595 !
596 ! IF GRIB2 read out post aviable fields xml file and post control file
597 !
598  if(grib == "grib2") then
599  btim=mpi_wtime()
600  call read_xml()
601  readxml_tim = readxml_tim + (mpi_wtime() - btim)
602  endif
603 !
604 ! LOOP OVER THE OUTPUT GRID(S). FIELD(S) AND OUTPUT GRID(S) ARE SPECIFIED
605 ! IN THE CONTROL FILE. WE PROCESS ONE GRID AND ITS FIELDS AT A TIME.
606 ! THAT'S WHAT THIS LOOP DOES.
607 !
608  icount_calmict = 0
609  first_grbtbl = .true.
610  npset = 0
611 !10 CONTINUE
612 !
613 ! READ CONTROL FILE DIRECTING WHICH FIELDS ON WHICH
614 ! LEVELS AND TO WHICH GRID TO INTERPOLATE DATA TO.
615 ! VARIABLE IEOF/=0 WHEN THERE ARE NO MORE GRIDS TO PROCESS.
616 !
617 ! -------- grib1 processing ---------------
618 ! ------------------
619 ! if (grib == "grib1") then !DO NOT REVERT TO GRIB1. GRIB1 NOT SUPPORTED ANYMORE
620 ! IEOF = 0
621 ! do while (ieof == 0)
622 ! CALL READCNTRL(kth,IEOF)
623 ! IF(ME == 0)THEN
624 ! WRITE(6,*)'POST: RETURN FROM READCNTRL. ', 'IEOF=',IEOF
625 ! ENDIF
626 !
627 ! PROCESS SELECTED FIELDS. FOR EACH SELECTED FIELD/LEVEL
628 ! WE GO THROUGH THE FOLLOWING STEPS:
629 ! (1) COMPUTE FIELD IF NEED BE
630 ! (2) WRITE FIELD TO OUTPUT FILE IN GRIB.
631 !
632 ! if (ieof == 0) then
633 ! CALL PROCESS(kth,kpv,th(1:kth),pv(1:kpv),iostatusD3D)
634 ! IF(ME == 0)THEN
635 ! WRITE(6,*)' '
636 ! WRITE(6,*)'WRFPOST: PREPARE TO PROCESS NEXT GRID'
637 ! ENDIF
638 ! endif
639 !
640 ! PROCESS NEXT GRID.
641 !
642 ! enddo
643 ! -------- grib2 processing ---------------
644 ! ------------------
645 ! elseif (grib == "grib2") then
646  if (me==0) write(0,*) ' in WRFPOST OUTFORM= ',grib
647  if (me==0) write(0,*) ' GRIB1 IS NOT SUPPORTED ANYMORE'
648  if (grib == "grib2") then
649  do while (npset < num_pset)
650  npset = npset+1
651  if (me==0) write(0,*)' in WRFPOST npset=',npset,' num_pset=',num_pset
652  CALL set_outflds(kth,th,kpv,pv)
653  if (me==0) write(0,*)' in WRFPOST size datapd',size(datapd)
654  if(allocated(datapd)) deallocate(datapd)
655  allocate(datapd(im,1:jend-jsta+1,nrecout+100))
656 !$omp parallel do private(i,j,k)
657  do k=1,nrecout+100
658  do j=1,jend+1-jsta
659  do i=1,im
660  datapd(i,j,k) = 0.
661  enddo
662  enddo
663  enddo
664  call get_postfilename(post_fname)
665  if (me==0) write(0,*)'post_fname=',trim(post_fname)
666  if (me==0) write(0,*)'get_postfilename,post_fname=',trim(post_fname), &
667  'npset=',npset, 'num_pset=',num_pset, &
668  'iSF_SURFACE_PHYSICS=',isf_surface_physics
669 !
670 ! PROCESS SELECTED FIELDS. FOR EACH SELECTED FIELD/LEVEL
671 ! WE GO THROUGH THE FOLLOWING STEPS:
672 ! (1) COMPUTE FIELD IF NEED BE
673 ! (2) WRITE FIELD TO OUTPUT FILE IN GRIB.
674 !
675  CALL process(kth,kpv,th(1:kth),pv(1:kpv),iostatusd3d)
676  IF(me == 0) WRITE(6,*)'WRFPOST: PREPARE TO PROCESS NEXT GRID'
677 !
678 ! write(0,*)'enter gribit2 before mpi_barrier'
679  call mpi_barrier(mpi_comm_comp,ierr)
680 
681 ! if(me==0)call w3tage('bf grb2 ')
682 ! write(0,*)'enter gribit2 after mpi barrier'
683  call gribit2(post_fname)
684  deallocate(datapd)
685  deallocate(fld_info)
686 !
687 ! PROCESS NEXT GRID.
688 !
689  enddo
690 
691  endif
692 !
693 !-------
694  call grib_info_finalize()
695 !
696  IF(me == 0) THEN
697  WRITE(6,*)' '
698  WRITE(6,*)'ALL GRIDS PROCESSED.'
699  WRITE(6,*)' '
700  ENDIF
701 !
702  call de_allocate
703 
704 ! GO TO 98
705  1000 CONTINUE
706 !exp call ext_ncd_ioclose ( DataHandle, Status )
707 !
708  IF(me == 0) THEN
709  print*, 'INITPOST_tim = ', initpost_tim
710  print*, 'MDLFLD_tim = ', etafld2_tim
711  print*, 'MDL2P_tim = ',eta2p_tim
712  print*, 'MDL2SIGMA_tim = ',mdl2sigma_tim
713  print*, 'MDL2AGL_tim = ',mdl2agl_tim
714  print*, 'SURFCE_tim = ',surfce2_tim
715  print*, 'CLDRAD_tim = ',cldrad_tim
716  print*, 'MISCLN_tim = ',miscln_tim
717  print*, 'MDL2STD_tim = ',mdl2std_tim
718  print*, 'FIXED_tim = ',fixed_tim
719  print*, 'MDL2THANDPV_tim = ',mdl2thandpv_tim
720  print*, 'CALRAD_WCLOUD_tim = ',calrad_wcloud_tim
721  print*, 'Total time = ',(mpi_wtime() - bbtim)
722  print*, 'Time for OUTPUT = ',time_output
723  print*, 'Time for READxml = ',readxml_tim
724  endif
725 !
726 ! END OF PROGRAM.
727 !
728 !
729 ! MPI_LAST WILL SHUTDOWN THE IO SERVER, IF IT EXISTS
730 !
731  CALL mpi_last
732 !
733 !
734  end if
735 !
736 !
737 !
738 ! call summary()
739  if (me == 0) CALL w3tage('UNIFIED_POST')
740  CALL mpi_finalize(ierr)
741 
742 
743  stop 0
744 
745  END
746