Update version info for release v4.6.1 (#2122)
[WRF.git] / external / io_grib1 / io_grib1.F
blob05864e26d0903ee41bdb4a500d114d19b9445d7d
1 !*-----------------------------------------------------------------------------
2 !*
3 !*  Todd Hutchinson
4 !*  WSI
5 !*  400 Minuteman Road
6 !*  Andover, MA     01810
7 !*  thutchinson@wsi.com
8 !*
9 !*-----------------------------------------------------------------------------
12 !* This io_grib1 API is designed to read WRF input and write WRF output data
13 !*   in grib version 1 format.  
17 module gr1_data_info
20 !* This module will hold data internal to this I/O implementation.  
21 !*   The variables will be accessible by all functions (provided they have a 
22 !*   "USE gr1_data_info" line).
25   integer                , parameter       :: FATAL            = 1
26   integer                , parameter       :: DEBUG            = 100
27   integer                , parameter       :: DateStrLen       = 19
29   integer                , parameter       :: firstFileHandle  = 8
30   integer                , parameter       :: maxFileHandles   = 30 
31   integer                , parameter       :: maxLevels        = 1000
32   integer                , parameter       :: maxSoilLevels    = 100
33   integer                , parameter       :: maxPressLevels   = 100
34   integer                , parameter       :: maxTurbLayers    = 100
35   integer                , parameter       :: maxDomains       = 500
37   logical ,      dimension(maxFileHandles) :: committed, opened, used
38   character*128, dimension(maxFileHandles) :: DataFile
39   integer,       dimension(maxFileHandles) :: FileFd
40   integer,       dimension(maxFileHandles) :: FileStatus
41   REAL,          dimension(maxLevels)      :: half_eta, full_eta
42   REAL,          dimension(maxSoilLevels)  :: soil_depth, soil_thickness
43   REAL,          dimension(maxPressLevels) :: press_levels
44   REAL,          dimension(maxTurbLayers)  :: turb_layer_bot
45   REAL,          dimension(maxTurbLayers)  :: turb_layer_top
46   character*24                             :: StartDate = ''
47   character*24                             :: InputProgramName = ''
48   integer                                  :: projection
49   integer                                  :: wg_grid_id
50   real                                     :: dx,dy
51   real                                     :: truelat1, truelat2
52   real                                     :: center_lat, center_lon
53   real                                     :: proj_central_lon
54   real                                     :: timestep
55   character,     dimension(:), pointer     :: grib_tables
56   logical                                  :: table_filled = .FALSE.
57   character,     dimension(:), pointer     :: grid_info
58   integer                                  :: full_xsize, full_ysize
59   integer, dimension(maxDomains)           :: domains = -1
60   integer                                  :: this_domain = 0
61   integer                                  :: max_domain = 0
62   
63   TYPE :: HandleVar
64      character, dimension(:), pointer      :: fileindex(:)
65      integer                               :: CurrentTime
66      integer                               :: NumberTimes
67      character (DateStrLen), dimension(:),pointer  :: Times(:)
68   ENDTYPE
69   TYPE (HandleVar), dimension(maxFileHandles) :: fileinfo
71   TYPE :: prevdata
72      integer :: fcst_secs_rainc
73      integer :: fcst_secs_rainnc
74      real, dimension(:,:), pointer         :: rainc, rainnc
75   END TYPE prevdata
76   TYPE (prevdata), DIMENSION(500)          :: lastdata
78   TYPE :: initdata
79      real,         dimension(:,:), pointer :: snod
80   END TYPE initdata
82   TYPE (initdata), dimension(maxDomains)   :: firstdata
84   TYPE :: prestype
85      real,         dimension(:,:,:), pointer :: vals
86      logical                                :: newtime
87      character*120                          :: lastDateStr
88   END TYPE prestype
90   character*120, dimension(maxDomains)     :: lastDateStr
92   TYPE (prestype), dimension(maxDomains)   :: pressure
93   TYPE (prestype), dimension(maxDomains)   :: geopotential
95   integer                                  :: center, subcenter, parmtbl
97   character(len=15000), dimension(firstFileHandle:maxFileHandles) :: td_output
98   character(len=15000), dimension(firstFileHandle:maxFileHandles) :: ti_output
100   logical                                  :: WrfIOnotInitialized = .true.
102 end module gr1_data_info
105 subroutine ext_gr1_ioinit(SysDepInfo,Status)
107   USE gr1_data_info
108   implicit none
109 #include "wrf_status_codes.h"
110 #include "wrf_io_flags.h"
111   CHARACTER*(*), INTENT(IN) :: SysDepInfo
112   integer ,intent(out) :: Status
113   integer :: i
114   integer :: size, istat
115   CHARACTER (LEN=300) :: wrf_err_message
117   call wrf_debug ( DEBUG , 'Entering ext_gr1_ioinit')
119   do i=firstFileHandle, maxFileHandles
120         used(i) = .false.
121         committed(i) = .false.
122         opened(i) = .false.
123         td_output(i) = ''
124         ti_output(i) = ''
125   enddo
126   domains(:) = -1
128   do i = 1, maxDomains
129     pressure(i)%newtime = .false.
130     pressure(i)%lastDateStr = ''
131     geopotential(i)%newtime = .false.
132     geopotential(i)%lastDateStr = ''
133     lastDateStr(i) = ''
134   enddo
136   lastdata%fcst_secs_rainc = 0
137   lastdata%fcst_secs_rainnc = 0
138   FileStatus(1:maxFileHandles) = WRF_FILE_NOT_OPENED
139   WrfIOnotInitialized = .false.
141   Status = WRF_NO_ERR
143   return
144 end subroutine ext_gr1_ioinit
146 !*****************************************************************************
148 subroutine ext_gr1_ioexit(Status)
150   USE gr1_data_info
151   implicit none
152 #include "wrf_status_codes.h"
153   integer istat
154   integer ,intent(out) :: Status
156   call wrf_debug ( DEBUG , 'Entering ext_gr1_ioexit')
158   if (table_filled) then
159      CALL free_gribmap(grib_tables)
160      DEALLOCATE(grib_tables, stat=istat)
161      table_filled = .FALSE.
162   endif
163   IF ( ASSOCIATED ( grid_info ) ) THEN
164     DEALLOCATE(grid_info, stat=istat)
165   ENDIF
166   NULLIFY(grid_info)
168   Status = WRF_NO_ERR
170   return
171 end subroutine ext_gr1_ioexit
173 !*****************************************************************************
175 SUBROUTINE ext_gr1_open_for_read_begin ( FileName , Comm_compute, Comm_io, &
176      SysDepInfo, DataHandle , Status )
178   USE gr1_data_info
179   IMPLICIT NONE
180 #include "wrf_status_codes.h"
181 #include "wrf_io_flags.h"
182   CHARACTER*(*) :: FileName
183   INTEGER ,       INTENT(IN)  :: Comm_compute , Comm_io
184   CHARACTER*(*) :: SysDepInfo
185   INTEGER ,       INTENT(OUT) :: DataHandle
186   INTEGER ,       INTENT(OUT) :: Status
187   integer                     :: ierr
188   integer                     :: size
189   integer                     :: idx
190   integer                     :: parmid
191   integer                     :: dpth_parmid
192   integer                     :: thk_parmid
193   integer                     :: leveltype
194   integer , DIMENSION(1000)   :: indices
195   integer                     :: numindices
196   real , DIMENSION(1000)      :: levels
197   real                        :: tmp
198   integer                     :: swapped
199   integer                     :: etaidx
200   integer                     :: grb_index
201   integer                     :: level1, level2
202   integer   :: tablenum
203   integer   :: stat
204   integer   :: endchar
205   integer   :: last_grb_index
206   CHARACTER (LEN=300) :: wrf_err_message
208   call wrf_debug ( DEBUG , 'Entering ext_gr1_open_for_read_begin')
210   CALL gr1_get_new_handle(DataHandle)
212   if (DataHandle .GT. 0) then
213      CALL open_file(TRIM(FileName), 'r', FileFd(DataHandle), ierr)
214      if (ierr .ne. 0) then
215         Status = WRF_ERR_FATAL_BAD_FILE_STATUS
216      else
217         opened(DataHandle) = .true.
218         DataFile(DataHandle) = TRIM(FileName)
219         FileStatus(DataHandle) = WRF_FILE_OPENED_NOT_COMMITTED
220      endif
221   else
222      Status = WRF_WARN_TOO_MANY_FILES
223      return
224   endif
226   ! Read the grib index file first
227   if (.NOT. table_filled) then
228      table_filled = .TRUE.
229      CALL GET_GRIB1_TABLES_SIZE(size)
230      ALLOCATE(grib_tables(1:size), STAT=ierr)
231      CALL LOAD_GRIB1_TABLES ("gribmap.txt", grib_tables, ierr)
232      if (ierr .ne. 0) then
233         DEALLOCATE(grib_tables)
234         WRITE( wrf_err_message , * ) &
235              'Could not open file gribmap.txt '
236         CALL wrf_error_fatal ( TRIM ( wrf_err_message ) )
237         Status = WRF_ERR_FATAL_BAD_FILE_STATUS
238         return
239      endif
240   endif
242   ! Begin by indexing file and reading metadata into structure.
243   CALL GET_FILEINDEX_SIZE(size)
244   ALLOCATE(fileinfo(DataHandle)%fileindex(1:size), STAT=ierr)
246   CALL ALLOC_INDEX_FILE(fileinfo(DataHandle)%fileindex(:))
247   CALL INDEX_FILE(FileFd(DataHandle),fileinfo(DataHandle)%fileindex(:))
249   ! Get times into Times variable
250   CALL GET_NUM_TIMES(fileinfo(DataHandle)%fileindex(:), &
251        fileinfo(DataHandle)%NumberTimes);
253   ALLOCATE(fileinfo(DataHandle)%Times(1:fileinfo(DataHandle)%NumberTimes), STAT=ierr)
254   do idx = 1,fileinfo(DataHandle)%NumberTimes
255      CALL GET_TIME(fileinfo(DataHandle)%fileindex(:),idx, &
256           fileinfo(DataHandle)%Times(idx))
257   enddo
259   ! CurrentTime starts as 0.  The first time in the file is 1.  So,
260   !   until set_time or get_next_time is called, the current time
261   !   is not set.
262   fileinfo(DataHandle)%CurrentTime = 0
264   CALL gr1_fill_eta_levels(fileinfo(DataHandle)%fileindex(:), &
265        FileFd(DataHandle), & 
266        grib_tables, "ZNW", full_eta)
267   CALL gr1_fill_eta_levels(fileinfo(DataHandle)%fileindex(:), FileFd(DataHandle), &
268        grib_tables, "ZNU", half_eta)
270   ! 
271   ! Now, get the soil levels
272   !
273   CALL GET_GRIB_PARAM(grib_tables, "ZS", center, subcenter, parmtbl, &
274        tablenum, dpth_parmid)
275   CALL GET_GRIB_PARAM(grib_tables,"DZS", center, subcenter, parmtbl, &
276        tablenum, thk_parmid)
277   if (dpth_parmid == -1) then
278      call wrf_message ('Error getting grib parameter')
279   endif
281   leveltype = 112
283   CALL GET_GRIB_INDICES(fileinfo(DataHandle)%fileindex(:),center, subcenter, parmtbl, &
284        dpth_parmid,"*",leveltype, &
285        -HUGE(1),-HUGE(1), -HUGE(1),-HUGE(1),indices,numindices)
287   last_grb_index = -1;
288   do idx = 1,numindices
289      CALL READ_GRIB(fileinfo(DataHandle)%fileindex(:), FileFd(DataHandle), &
290           indices(idx), soil_depth(idx))
291      !
292      ! Now read the soil thickenesses
293      !
294      CALL GET_LEVEL1(fileinfo(DataHandle)%fileindex(:),indices(idx),level1)
295      CALL GET_LEVEL2(fileinfo(DataHandle)%fileindex(:),indices(idx),level2)
296      CALL GET_GRIB_INDEX_GUESS(fileinfo(DataHandle)%fileindex(:), &
297           center, subcenter, parmtbl, thk_parmid,"*",leveltype, &
298           level1,level2,-HUGE(1),-HUGE(1), last_grb_index+1, grb_index)
299      CALL READ_GRIB(fileinfo(DataHandle)%fileindex(:),FileFd(DataHandle),grb_index, &
300           soil_thickness(idx))
302      last_grb_index = grb_index
303   enddo
304   
307   !
308   ! Fill up any variables that need to be retrieved from Metadata
309   !
310   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), 'PROGRAM_NAME', "none", &
311        "none", InputProgramName, stat)
312   if (stat /= 0) then
313      CALL wrf_debug (DEBUG , "PROGRAM_NAME not found in input METADATA")
314   else 
315      endchar = SCAN(InputProgramName," ")
316      InputProgramName = InputProgramName(1:endchar)
317   endif
319   call wrf_debug ( DEBUG , 'Exiting ext_gr1_open_for_read_begin')
321   RETURN
322 END SUBROUTINE ext_gr1_open_for_read_begin
324 !*****************************************************************************
326 SUBROUTINE ext_gr1_open_for_read_commit( DataHandle , Status )
328   USE gr1_data_info
329   IMPLICIT NONE
330 #include "wrf_status_codes.h"
331 #include "wrf_io_flags.h"
332   character(len=1000) :: msg
333   INTEGER ,       INTENT(IN ) :: DataHandle
334   INTEGER ,       INTENT(OUT) :: Status
336   call wrf_debug ( DEBUG , 'Entering ext_gr1_open_for_read_commit')
338   Status = WRF_NO_ERR
339   if(WrfIOnotInitialized) then
340     Status = WRF_IO_NOT_INITIALIZED
341     write(msg,*) 'ext_gr1_ioinit was not called ',__FILE__,', line', __LINE__
342     call wrf_debug ( FATAL , msg)
343     return
344   endif
345   committed(DataHandle) = .true.
346   FileStatus(DataHandle) = WRF_FILE_OPENED_FOR_READ
348   Status = WRF_NO_ERR
350   RETURN
351 END SUBROUTINE ext_gr1_open_for_read_commit
353 !*****************************************************************************
355 SUBROUTINE ext_gr1_open_for_read ( FileName , Comm_compute, Comm_io, &
356      SysDepInfo, DataHandle , Status )
358   USE gr1_data_info
359   IMPLICIT NONE
360 #include "wrf_status_codes.h"
361 #include "wrf_io_flags.h"
362   CHARACTER*(*) :: FileName
363   INTEGER ,       INTENT(IN)  :: Comm_compute , Comm_io
364   CHARACTER*(*) :: SysDepInfo
365   INTEGER ,       INTENT(OUT) :: DataHandle
366   INTEGER ,       INTENT(OUT) :: Status
369   call wrf_debug ( DEBUG , 'Entering ext_gr1_open_for_read')
371   DataHandle = 0   ! dummy setting to quiet warning message
372   CALL ext_gr1_open_for_read_begin( FileName, Comm_compute, Comm_io, &
373        SysDepInfo, DataHandle, Status )
374   IF ( Status .EQ. WRF_NO_ERR ) THEN
375      FileStatus(DataHandle) = WRF_FILE_OPENED_NOT_COMMITTED
376      CALL ext_gr1_open_for_read_commit( DataHandle, Status )
377   ENDIF
378   return
380   RETURN  
381 END SUBROUTINE ext_gr1_open_for_read
383 !*****************************************************************************
385 SUBROUTINE ext_gr1_open_for_write_begin(FileName, Comm, IOComm, SysDepInfo, &
386      DataHandle, Status)
387   
388   USE gr1_data_info
389   implicit none
390 #include "wrf_status_codes.h"
391 #include "wrf_io_flags.h"
393   character*(*)        ,intent(in)  :: FileName
394   integer              ,intent(in)  :: Comm
395   integer              ,intent(in)  :: IOComm
396   character*(*)        ,intent(in)  :: SysDepInfo
397   integer              ,intent(out) :: DataHandle
398   integer              ,intent(out) :: Status
399   integer :: ierr
400   CHARACTER (LEN=300) :: wrf_err_message
401   integer             :: size
403   call wrf_debug ( DEBUG , 'Entering ext_gr1_open_for_write_begin')
405   if (.NOT. table_filled) then
406      table_filled = .TRUE.
407      CALL GET_GRIB1_TABLES_SIZE(size)
408      ALLOCATE(grib_tables(1:size), STAT=ierr)
409      CALL LOAD_GRIB1_TABLES ("gribmap.txt", grib_tables, ierr)
410      if (ierr .ne. 0) then
411         DEALLOCATE(grib_tables)
412         WRITE( wrf_err_message , * ) &
413              'Could not open file gribmap.txt '
414         CALL wrf_error_fatal ( TRIM ( wrf_err_message ) )
415         Status = WRF_ERR_FATAL_BAD_FILE_STATUS
416         return
417      endif
418   endif
420   Status = WRF_NO_ERR
421   CALL gr1_get_new_handle(DataHandle)
422   if (DataHandle .GT. 0) then
423      CALL open_file(TRIM(FileName), 'w', FileFd(DataHandle), ierr)
424      if (ierr .ne. 0) then
425         Status = WRF_WARN_WRITE_RONLY_FILE
426      else
427         opened(DataHandle) = .true.
428         DataFile(DataHandle) = TRIM(FileName)
429         FileStatus(DataHandle) = WRF_FILE_OPENED_NOT_COMMITTED
430      endif
431      committed(DataHandle) = .false.
432      td_output(DataHandle) = ''
433   else
434      Status = WRF_WARN_TOO_MANY_FILES
435   endif
437   RETURN  
438 END SUBROUTINE ext_gr1_open_for_write_begin
440 !*****************************************************************************
442 SUBROUTINE ext_gr1_open_for_write_commit( DataHandle , Status )
444   USE gr1_data_info
445   IMPLICIT NONE
446 #include "wrf_status_codes.h"
447 #include "wrf_io_flags.h"
448   INTEGER ,       INTENT(IN ) :: DataHandle
449   INTEGER ,       INTENT(OUT) :: Status
451   call wrf_debug ( DEBUG , 'Entering ext_gr1_open_for_write_commit')
453   IF ( opened( DataHandle ) ) THEN
454     IF ( used( DataHandle ) ) THEN
455       committed(DataHandle) = .true.
456       FileStatus(DataHandle) = WRF_FILE_OPENED_FOR_WRITE
457     ENDIF
458   ENDIF
460   Status = WRF_NO_ERR
462   RETURN  
463 END SUBROUTINE ext_gr1_open_for_write_commit
465 !*****************************************************************************
467 subroutine ext_gr1_inquiry (Inquiry, Result, Status)
468   use gr1_data_info
469   implicit none
470 #include "wrf_status_codes.h"
471   character *(*), INTENT(IN)    :: Inquiry
472   character *(*), INTENT(OUT)   :: Result
473   integer        ,INTENT(INOUT) :: Status
474   SELECT CASE (Inquiry)
475   CASE ("RANDOM_WRITE","RANDOM_READ")
476      Result='ALLOW'
477   CASE ("SEQUENTIAL_WRITE","SEQUENTIAL_READ")
478      Result='NO'
479   CASE ("OPEN_READ", "OPEN_WRITE", "OPEN_COMMIT_WRITE")
480      Result='REQUIRE'
481   CASE ("OPEN_COMMIT_READ","PARALLEL_IO")
482      Result='NO'
483   CASE ("SELF_DESCRIBING","SUPPORT_METADATA","SUPPORT_3D_FIELDS")
484      Result='YES'
485   CASE ("MEDIUM")
486      Result ='FILE'
487   CASE DEFAULT
488      Result = 'No Result for that inquiry!'
489   END SELECT
490   Status=WRF_NO_ERR
491   return
492 end subroutine ext_gr1_inquiry
494 !*****************************************************************************
496 SUBROUTINE ext_gr1_inquire_opened ( DataHandle, FileName , FileStat, Status )
498   USE gr1_data_info
499   IMPLICIT NONE
500 #include "wrf_status_codes.h"
501 #include "wrf_io_flags.h"
502   INTEGER ,       INTENT(IN)  :: DataHandle
503   CHARACTER*(*) :: FileName
504   INTEGER ,       INTENT(OUT) :: FileStat
505   INTEGER ,       INTENT(OUT) :: Status
507   call wrf_debug ( DEBUG , 'Entering ext_gr1_inquire_opened')
509   FileStat = WRF_NO_ERR
510   if ((DataHandle .ge. firstFileHandle) .and. &
511        (DataHandle .le. maxFileHandles)) then
512      FileStat = FileStatus(DataHandle)
513   else
514      FileStat = WRF_FILE_NOT_OPENED
515   endif
516   
517   Status = FileStat
519   RETURN
520 END SUBROUTINE ext_gr1_inquire_opened
522 !*****************************************************************************
524 SUBROUTINE ext_gr1_ioclose ( DataHandle, Status )
526   USE gr1_data_info
527   IMPLICIT NONE
528 #include "wrf_status_codes.h"
529   INTEGER DataHandle, Status
530   INTEGER istat
531   INTEGER ierr
532   character(len=1000) :: outstring
533   character :: lf
534   lf=char(10)
535      
536   call wrf_debug ( DEBUG , 'Entering ext_gr1_ioclose')
538   Status = WRF_NO_ERR
540   CALL write_file(FileFd(DataHandle), lf//'<METADATA>'//lf,ierr)
541   outstring = &
542        '<!-- The following are fields that were supplied to the WRF I/O API.'//lf//&
543        'Many variables (but not all) are redundant with the variables within '//lf//&
544        'the grib headers.  They are stored here, as METADATA, so that the '//lf//&
545        'WRF I/O API has simple access to these variables.-->'
546   CALL write_file(FileFd(DataHandle), trim(outstring), ierr)
547   if (trim(ti_output(DataHandle)) /= '') then
548      CALL write_file(FileFd(DataHandle), trim(ti_output(DataHandle)), ierr)
549      CALL write_file(FileFd(DataHandle), lf, ierr)
550   endif
551   if (trim(td_output(DataHandle)) /= '') then
552      CALL write_file(FileFd(DataHandle), trim(td_output(DataHandle)), ierr)
553      CALL write_file(FileFd(DataHandle), lf, ierr)
554   endif
555   CALL write_file(FileFd(DataHandle), '</METADATA>'//lf,ierr)
556   ti_output(DataHandle) = ''
557   td_output(DataHandle) = ''
558   if (ierr .ne. 0) then
559      Status = WRF_WARN_WRITE_RONLY_FILE
560   endif
561   CALL close_file(FileFd(DataHandle))
563   used(DataHandle) = .false.
565   RETURN
566 END SUBROUTINE ext_gr1_ioclose
568 !*****************************************************************************
570 SUBROUTINE ext_gr1_write_field( DataHandle , DateStrIn , VarName , &
571      Field , FieldType , Comm , IOComm, &
572      DomainDesc , MemoryOrder , Stagger , &
573      DimNames , &
574      DomainStart , DomainEnd , &
575      MemoryStart , MemoryEnd , &
576      PatchStart , PatchEnd , &
577      Status )
579   USE gr1_data_info
580   IMPLICIT NONE
581 #include "wrf_status_codes.h"
582 #include "wrf_io_flags.h"
583 #include "wrf_projection.h"
584   INTEGER ,       INTENT(IN)    :: DataHandle 
585   CHARACTER*(*) :: DateStrIn
586   CHARACTER(DateStrLen) :: DateStr
587   CHARACTER*(*) :: VarName
588   CHARACTER*120 :: OutName
589   CHARACTER(120) :: TmpVarName
590   integer                       ,intent(in)    :: FieldType
591   integer                       ,intent(inout) :: Comm
592   integer                       ,intent(inout) :: IOComm
593   integer                       ,intent(in)    :: DomainDesc
594   character*(*)                 ,intent(in)    :: MemoryOrder
595   character*(*)                 ,intent(in)    :: Stagger
596   character*(*) , dimension (*) ,intent(in)    :: DimNames
597   integer ,dimension(*)         ,intent(in)    :: DomainStart, DomainEnd
598   integer ,dimension(*)         ,intent(in)    :: MemoryStart, MemoryEnd
599   integer ,dimension(*)         ,intent(in)    :: PatchStart,  PatchEnd
600   integer                       ,intent(out)   :: Status
601   integer                                      :: ierror
602   character (120)                         :: msg
603   integer :: xsize, ysize, zsize
604   integer :: x, y, z
605   integer :: x_start,x_end,y_start,y_end,z_start,z_end,ndim
606   integer :: idx
607   integer :: proj_center_flag
608   logical :: vert_stag = .false.
609   integer :: levelnum
610   real, DIMENSION(:,:), POINTER :: data,tmpdata
611   integer, DIMENSION(:), POINTER :: mold
612   integer :: istat
613   integer :: accum_period
614   integer :: size
615   integer, dimension(1000) :: level1, level2
616   real, DIMENSION( 1:1,MemoryStart(1):MemoryEnd(1), &
617                    MemoryStart(2):MemoryEnd(2), &
618                    MemoryStart(3):MemoryEnd(3) ) :: Field
619   real    :: fcst_secs
620   logical :: soil_layers, fraction, is_press_levels,is_turb_layers
621   integer :: vert_unit
622   integer :: abc(2,2,2)
623   integer :: def(8)
624   logical :: output = .true.
625   integer :: idx1, idx2, idx3
626   logical :: new_domain
627   real    :: region_center_lat, region_center_lon
628   integer :: dom_xsize, dom_ysize;
629   integer :: ierr
630   logical :: already_have_domain
632   call wrf_debug ( DEBUG , 'Entering ext_gr1_write_field for parameter'//VarName)
634   !
635   ! If DateStr is all 0's, we reset it to StartDate (if StartDate exists).  
636   !   For some reason, 
637   !   in idealized simulations, StartDate is 0001-01-01_00:00:00 while
638   !   the first DateStr is 0000-00-00_00:00:00.  
639   !
640   if (DateStrIn .eq. '0000-00-00_00:00:00') then
641      if (StartDate .ne. '') then
642         DateStr = TRIM(StartDate)
643      else
644         DateStr = '0001-01-01_00:00:00'
645      endif
646   else
647      DateStr = DateStrIn
648   endif
650   !
651   ! Check if this is a domain that we haven't seen yet.  If so, add it to 
652   !   the list of domains.
653   !
654   new_domain = .false.
655   already_have_domain = .false.
656   do idx = 1, max_domain
657      if (this_domain .eq. domains(idx)) then
658         already_have_domain = .true.
659      endif
660   enddo
661   if (.NOT. already_have_domain) then
662      max_domain = max_domain + 1
663      domains(max_domain) = this_domain
664      new_domain = .true.
665   endif
667   !
668   ! If the time has changed, we open a new file.  This is a kludge to get
669   !   around slowness in WRF that occurs when opening a new data file the
670   !   standard way.
671   !
672 #ifdef GRIB_ONE_TIME_PER_FILE
673   if (lastDateStr(this_domain) .ne. DateStr) then
674      write(DataFile(DataHandle),'(A8,i2.2,A1,A19)') 'wrfout_d',this_domain,'_',DateStr
675      call ext_gr1_ioclose ( DataHandle, Status )
676      CALL open_file(TRIM(DataFile(DataHandle)), 'w', FileFd(DataHandle), ierr)
677      if (ierr .ne. 0) then
678         print *,'Could not open new file: ',DataFile(DataHandle)
679         print *,'  Appending to old file.'
680      else
681         ! Just set used back to .true. here, since ioclose set it to false.
682         used(DataHandle) = .true.
683      endif
684      td_output(DataHandle) = ''
685   endif
686   lastDateStr(this_domain) = DateStr
687 #endif
689   output = .true.
690   zsize = 1
691   xsize = 1
692   ysize = 1
693   OutName = VarName
694   is_press_levels = .false.
695   is_turb_layers = .false.
696   soil_layers = .false.
697   fraction = .false.
699   ! First, handle then special cases for the boundary data.
701   CALL get_dims(MemoryOrder, PatchStart, PatchEnd, ndim, x_start, x_end, &
702        y_start, y_end,z_start,z_end)
703   xsize = x_end - x_start + 1
704   ysize = y_end - y_start + 1
705   zsize = z_end - z_start + 1
707   do idx = 1, len(MemoryOrder)
708      if ((MemoryOrder(idx:idx) .eq. 'Z') .and. &
709           (DimNames(idx) .eq. 'soil_layers_stag')) then
710         soil_layers = .true.
711      else if ((MemoryOrder(idx:idx) .eq. 'Z') .and. &
712           (DimNames(idx) .eq. 'num_press_levels_stag')) then
713         is_press_levels = .true.
714      else if ((MemoryOrder(idx:idx) .eq. 'Z') .and. &
715           (DimNames(idx) .eq. 'num_turb_layers')) then
716         is_turb_layers = .true.
717      else if ((OutName .eq. 'LANDUSEF') .or. (OutName .eq. 'SOILCBOT') .or. &
718           (OutName .eq. 'SOILCTOP')) then
719         fraction = .true.
720      endif
721   enddo
723   if (.not. ASSOCIATED(grid_info)) then
724      CALL get_grid_info_size(size)
725      ALLOCATE(grid_info(1:size), STAT=istat)
726      if (istat .eq. -1) then
727         DEALLOCATE(grid_info)
728         Status = WRF_ERR_FATAL_BAD_FILE_STATUS
729         return
730      endif
731   endif
732      
734   if (new_domain) then
735      ALLOCATE(firstdata(this_domain)%snod(xsize,ysize))
736      firstdata(this_domain)%snod(:,:) = 0.0
737      ALLOCATE(lastdata(this_domain)%rainc(xsize,ysize))
738      lastdata(this_domain)%rainc(:,:) = 0.0
739      ALLOCATE(lastdata(this_domain)%rainnc(xsize,ysize))
740      lastdata(this_domain)%rainnc(:,:) = 0.0
741   endif
743   if (zsize .eq. 0) then 
744      zsize = 1
745   endif
747   ALLOCATE(data(1:xsize,1:ysize), STAT=istat)
748   ALLOCATE(mold(1:ysize), STAT=istat)
749   ALLOCATE(tmpdata(1:xsize,1:ysize), STAT=istat)
751   if (OutName .eq. 'ZNU') then
752      do idx = 1, zsize
753         half_eta(idx) = Field(1,idx,1,1)
754      enddo
755   endif
757   if (OutName .eq. 'ZNW') then
758      do idx = 1, zsize
759         full_eta(idx) = Field(1,idx,1,1)
760      enddo
761   endif
763   if (OutName .eq. 'ZS') then
764      do idx = 1, zsize
765         soil_depth(idx) = Field(1,idx,1,1)
766      enddo
767   endif
769   if (OutName .eq. 'DZS') then
770      do idx = 1, zsize
771         soil_thickness(idx) = Field(1,idx,1,1)
772      enddo
773   endif
775   if (OutName .eq. 'P_PL') then
776      do idx = 1, zsize
777         press_levels(idx) = Field(1,idx,1,1)
778      enddo
779   endif
781   if (OutName .eq. 'AFWA_TLYRBOT') then
782      do idx = 1, zsize
783         turb_layer_bot(idx) = Field(1,idx,1,1)
784      enddo
785   endif
787   if (OutName .eq. 'AFWA_TLYRTOP') then
788      do idx = 1, zsize
789         turb_layer_top(idx) = Field(1,idx,1,1)
790      enddo
791   endif
792   if ((xsize .lt. 1) .or. (ysize .lt. 1)) then
793      write(msg,*) 'Cannot output field with memory order: ', &
794           MemoryOrder,Varname
795      call wrf_message(msg)
796      return
797   endif
798      
799   call get_vert_stag(OutName,Stagger,vert_stag)
801   do idx = 1, zsize
802      call gr1_get_levels(OutName, idx, zsize, soil_layers, vert_stag, fraction, &
803           is_press_levels, is_turb_layers, vert_unit, level1(idx), level2(idx))
804   enddo
806   ! 
807   ! Get the center lat/lon for the area being output.  For some cases (such
808   !    as for boundary areas, the center of the area is different from the
809   !    center of the model grid.
810   !
811   if (index(Stagger,'X') .le. 0) then
812      dom_xsize = full_xsize - 1
813   else
814      dom_xsize = full_xsize
815   endif
816   if (index(Stagger,'Y') .le. 0) then
817      dom_ysize = full_ysize - 1
818   else
819      dom_ysize = full_ysize
820   endif
822   !
823   ! Handle case of polare stereographic centered on pole.  In that case,
824   !   always set center lon to be the projection central longitude.
825   !
826   if ((projection .eq. WRF_POLAR_STEREO) .AND. &
827        (abs(center_lat - 90.0) < 0.01)) then
828      center_lon = proj_central_lon
829   endif
831   CALL get_region_center(MemoryOrder, projection, center_lat, center_lon, &
832        dom_xsize, dom_ysize, dx, dy, proj_central_lon, proj_center_flag, &
833        truelat1, truelat2, xsize, ysize, region_center_lat, region_center_lon)
835   if ( .not. opened(DataHandle)) then
836      Status = WRF_WARN_FILE_NOT_OPENED
837      return
838   endif
841   if (opened(DataHandle) .and. committed(DataHandle)) then
844 #ifdef OUTPUT_FULL_PRESSURE
846      ! 
847      ! The following is a kludge to output full pressure instead of the two 
848      !  fields of base-state pressure and pressure perturbation.
849      !
850      ! This code can be turned on by adding -DOUTPUT_FULL_PRESSURE to the 
851      !  compile line
852      !
853      
854      if ((OutName .eq. 'P') .or. (OutName.eq.'PB')) then
855         do idx = 1, len(MemoryOrder)
856             if (MemoryOrder(idx:idx) .eq. 'X') then
857               idx1=idx
858            endif
859            if (MemoryOrder(idx:idx) .eq. 'Y') then
860               idx2=idx
861            endif
862            if (MemoryOrder(idx:idx) .eq. 'Z') then
863               idx3=idx
864            endif
865         enddo
867         ! 
868         ! Allocate space for pressure values (this variable holds 
869         !   base-state pressure or pressure perturbation to be used 
870         !   later to sum base-state and perturbation pressure to get full 
871         !   pressure).
872         !
874         if (.not. ASSOCIATED(pressure(this_domain)%vals)) then
875            ALLOCATE(pressure(this_domain)%vals(MemoryStart(1):MemoryEnd(1), &
876                 MemoryStart(2):MemoryEnd(2),MemoryStart(3):MemoryEnd(3)))
877         endif
878         if (DateStr .NE. &
879              pressure(this_domain)%lastDateStr) then
880            pressure(this_domain)%newtime = .true.
881         endif
882         if (pressure(this_domain)%newtime) then
883            pressure(this_domain)%vals = Field(1,:,:,:)
884            pressure(this_domain)%newtime = .false.
885            output = .false.
886         else 
887            output = .true.
888         endif
889         pressure(this_domain)%lastDateStr=DateStr
890      endif
891 #endif
893 #ifdef OUTPUT_FULL_GEOPOTENTIAL
895      ! 
896      ! The following is a kludge to output full geopotential height instead 
897      !  of the two fields of base-state geopotential and perturbation 
898      !  geopotential.
899      !
900      ! This code can be turned on by adding -DOUTPUT_FULL_GEOPOTENTIAL to the 
901      !  compile line
902      !
903      
904      if ((OutName .eq. 'PHB') .or. (OutName.eq.'PH')) then
905         do idx = 1, len(MemoryOrder)
906             if (MemoryOrder(idx:idx) .eq. 'X') then
907               idx1=idx
908            endif
909            if (MemoryOrder(idx:idx) .eq. 'Y') then
910               idx2=idx
911            endif
912            if (MemoryOrder(idx:idx) .eq. 'Z') then
913               idx3=idx
914            endif
915         enddo
917         ! 
918         ! Allocate space for geopotential values (this variable holds 
919         !   geopotential to be used 
920         !   later to sum base-state and perturbation to get full 
921         !   geopotential).
922         !
924         if (.not. ASSOCIATED(geopotential(this_domain)%vals)) then
925            ALLOCATE(geopotential(this_domain)%vals(MemoryStart(1):MemoryEnd(1), &
926                 MemoryStart(2):MemoryEnd(2),MemoryStart(3):MemoryEnd(3)))
927         endif
928         if (DateStr .NE. &
929              geopotential(this_domain)%lastDateStr) then
930            geopotential(this_domain)%newtime = .true.
931         endif
932         if (geopotential(this_domain)%newtime) then
933            geopotential(this_domain)%vals = Field(1,:,:,:)
934            geopotential(this_domain)%newtime = .false.
935            output = .false.
936         else 
937            output = .true.
938         endif
939         geopotential(this_domain)%lastDateStr=DateStr
940      endif
941 #endif
943      if (output) then 
944         if (StartDate == '') then
945            StartDate = DateStr
946         endif
947         CALL geth_idts(DateStr,StartDate,fcst_secs)
948         
949         if (center_lat .lt. 0) then
950            proj_center_flag = 2
951         else
952            proj_center_flag = 1
953         endif
954          
955         do z = 1, zsize
956            SELECT CASE (MemoryOrder)
957            CASE ('XYZ')
958               data = Field(1,1:xsize,1:ysize,z)
959            CASE ('XZY')
960               data = Field(1,1:xsize,z,1:ysize)
961            CASE ('YXZ')
962               do x = 1,xsize
963                  do y = 1,ysize
964                     data(x,y) = Field(1,y,x,z)
965                  enddo
966               enddo
967            CASE ('YZX')
968               do x = 1,xsize
969                  do y = 1,ysize
970                     data(x,y) = Field(1,y,z,x)
971                  enddo
972               enddo
973            CASE ('ZXY')
974               data = Field(1,z,1:xsize,1:ysize)
975            CASE ('ZYX')
976               do x = 1,xsize
977                  do y = 1,ysize
978                     data(x,y) = Field(1,z,y,x)
979                  enddo
980               enddo
981            CASE ('XY')
982               data = Field(1,1:xsize,1:ysize,1)
983            CASE ('YX')
984               do x = 1,xsize
985                  do y = 1,ysize
986                     data(x,y) = Field(1,y,x,1)
987                  enddo
988               enddo
990            CASE ('XSZ')
991               do x = 1,xsize
992                  do y = 1,ysize
993                     data(x,y) = Field(1,y,z,x)
994                  enddo
995               enddo
996            CASE ('XEZ')
997               do x = 1,xsize
998                  do y = 1,ysize
999                     data(x,y) = Field(1,y,z,x)
1000                  enddo
1001               enddo
1002            CASE ('YSZ')
1003               do x = 1,xsize
1004                  do y = 1,ysize
1005                     data(x,y) = Field(1,x,z,y)
1006                  enddo
1007               enddo
1008            CASE ('YEZ')
1009               do x = 1,xsize
1010                  do y = 1,ysize
1011                     data(x,y) = Field(1,x,z,y)
1012                  enddo
1013               enddo
1015            CASE ('XS')
1016               do x = 1,xsize
1017                  do y = 1,ysize
1018                     data(x,y) = Field(1,y,x,1)
1019                  enddo
1020               enddo
1021            CASE ('XE')
1022               do x = 1,xsize
1023                  do y = 1,ysize
1024                     data(x,y) = Field(1,y,x,1)
1025                  enddo
1026               enddo
1027            CASE ('YS')
1028               do x = 1,xsize
1029                  do y = 1,ysize
1030                     data(x,y) = Field(1,x,y,1)
1031                  enddo
1032               enddo
1033            CASE ('YE')
1034               do x = 1,xsize
1035                  do y = 1,ysize
1036                     data(x,y) = Field(1,x,y,1)
1037                  enddo
1038               enddo
1040            CASE ('Z')
1041               data(1,1) = Field(1,z,1,1)
1042            CASE ('z')
1043               data(1,1) = Field(1,z,1,1)
1044            CASE ('C')
1045               data = Field(1,1:xsize,1:ysize,z)
1046            CASE ('c')
1047               data = Field(1,1:xsize,1:ysize,z)
1048            CASE ('0')
1049               data(1,1) = Field(1,1,1,1)
1050            END SELECT
1052            ! 
1053            ! Here, we convert any integer fields to real
1054            !
1055            if (FieldType == WRF_INTEGER) then
1056               mold = 0
1057               do idx=1,xsize
1058                  !
1059                  ! The parentheses around data(idx,:) are needed in order
1060                  !   to fix a bug with transfer with the xlf compiler on NCAR's
1061                  !   IBM (bluesky).
1062                  !
1063                  data(idx,:)=transfer((data(idx,:)),mold)
1064               enddo
1065            endif
1066            ! 
1067            ! Here, we do any necessary conversions to the data.
1068            !
1069            
1070            ! Potential temperature is sometimes passed in as perturbation 
1071            !   potential temperature (i.e., POT-300).
1072            !
1073            if (OutName == 'T') then
1074               data = data + 300
1075            endif
1077            ! 
1078            ! For precip, we setup the accumulation period, and output a precip
1079            !    rate for time-step precip.
1080            !
1081            if (OutName .eq. 'RAINNCV') then
1082               ! Convert time-step precip to precip rate.
1083               data = data/timestep
1084               accum_period = 0
1085            else
1086               accum_period = 0
1087            endif
1089            if ((OutName .eq. 'AFWA_TOTPRECIP') .or. &
1090                (OutName .eq. 'AFWA_RAIN') .or. &
1091                (OutName .eq. 'AFWA_FZRA') .or. &
1092                (OutName .eq. 'AFWA_SNOW') .or. &
1093                (OutName .eq. 'AFWA_SNOWFALL') .or. &
1094                (OutName .eq. 'RAINC') .or. &
1095                (OutName .eq. 'AFWA_ICE')) then
1096                accum_period = fcst_secs
1097            endif
1099 #ifdef OUTPUT_FULL_PRESSURE
1100            !
1101            ! Computation of full-pressure off by default since there are 
1102            !  uses for base-state and perturbation (i.e., restarts
1103            !
1104             if ((OutName .eq. 'P') .or. (OutName.eq.'PB')) then
1105                if (idx3 .eq. 1) then
1106                   data = data + &
1107                        pressure(this_domain)%vals(z, &
1108                        patchstart(2):patchend(2),patchstart(3):patchend(3))
1109                elseif (idx3 .eq. 2) then
1110                   data = data + &
1111                        pressure(this_domain)%vals(patchstart(1):patchend(1), &
1112                        z,patchstart(3):patchend(3))
1113                elseif (idx3 .eq. 3) then
1114                   data = data + &
1115                        pressure(this_domain)%vals(patchstart(1):patchend(1), &
1116                        patchstart(2):patchend(2),z)
1117                else
1118                   call wrf_message ('error in idx3, continuing')
1119                endif
1121                OutName = 'P'
1122             endif
1123 #endif
1125 #ifdef OUTPUT_FULL_GEOPOTENTIAL
1126            !
1127            ! Computation of full-geopotential off by default since there are 
1128            !  uses for base-state and perturbation (i.e., restarts
1129            !
1130             if ((OutName .eq. 'PHB') .or. (OutName.eq.'PH')) then
1131                if (idx3 .eq. 1) then
1132                   data = data + &
1133                        geopotential(this_domain)%vals(z, &
1134                        patchstart(2):patchend(2),patchstart(3):patchend(3))
1135                elseif (idx3 .eq. 2) then
1136                   data = data + &
1137                        geopotential(this_domain)%vals(patchstart(1):patchend(1), &
1138                        z,patchstart(3):patchend(3))
1139                elseif (idx3 .eq. 3) then
1140                   data = data + &
1141                        geopotential(this_domain)%vals(patchstart(1):patchend(1), &
1142                        patchstart(2):patchend(2),z)
1143                else
1144                   call wrf_message ('error in idx3, continuing')
1145                endif
1147                OutName = 'PHP'
1148             endif
1149 #endif
1151            !
1152            !    Output current level
1153            !
1154            CALL load_grid_info(OutName, StartDate, vert_unit, level1(z), &
1155                 level2(z), fcst_secs, accum_period, wg_grid_id, projection, &
1156                 xsize, ysize, region_center_lat, region_center_lon, dx, dy, &
1157                 proj_central_lon, proj_center_flag, truelat1, truelat2, &
1158                 grib_tables, grid_info)
1159            
1160            !
1161            ! Here, we copy data to a temporary array.  After write_grib,
1162            !    we copy back from the temporary array to the permanent
1163            !    array.  write_grib modifies data.  For certain fields that
1164            !    we use below, we want the original (unmodified) data 
1165            !    values.  This kludge assures that we have the original
1166            !    values.
1167            !
1169            if ((OutName .eq. 'RAINC') .or. (OutName .eq. 'RAINNC') .or. &
1170                 (OutName .eq. 'SNOWH')) then
1171               tmpdata(:,:) = data(:,:)
1172            endif
1174            CALL write_grib(grid_info, FileFd(DataHandle), data)
1176            if ((OutName .eq. 'RAINC') .or. (OutName .eq. 'RAINNC') .or. &
1177                 (OutName .eq. 'SNOWH')) then
1178               data(:,:) = tmpdata(:,:)
1179            endif
1181            CALL free_grid_info(grid_info)
1182            
1183            !
1184            ! If this is the total accumulated rain, call write_grib again 
1185            !   to output the accumulation since the last output time as well.
1186            !   This is somewhat of a kludge to meet the requirements of PF.
1187            !
1188            if ((OutName .eq. 'RAINC') .or. (OutName .eq. 'RAINNC') .or. &
1189                 (OutName .eq. 'SNOWH')) then
1190               if (OutName .eq. 'RAINC') then
1191                  data(:,:) = data(:,:) - lastdata(this_domain)%rainc(:,:)
1192                  lastdata(this_domain)%rainc(:,:) = tmpdata(:,:)
1193                  accum_period = fcst_secs - &
1194                       lastdata(this_domain)%fcst_secs_rainc
1195                  lastdata(this_domain)%fcst_secs_rainc = fcst_secs
1196                  TmpVarName = 'ACPCP'
1197               else if (OutName .eq. 'RAINNC') then
1198                  tmpdata(:,:) = data(:,:)
1199                  data(:,:) = data(:,:) - lastdata(this_domain)%rainnc(:,:)
1200                  lastdata(this_domain)%rainnc(:,:) = tmpdata(:,:)
1201                  accum_period = fcst_secs - &
1202                       lastdata(this_domain)%fcst_secs_rainnc
1203                  lastdata(this_domain)%fcst_secs_rainnc = fcst_secs
1204                  TmpVarName = 'NCPCP'
1205               else if (OutName .eq. 'SNOWH') then
1206                  if (fcst_secs .eq. 0) then
1207                     firstdata(this_domain)%snod(:,:) = data(:,:)
1208                  endif
1209                  data(:,:) = data(:,:) - firstdata(this_domain)%snod(:,:)
1210                  TmpVarName = 'SNOWCU'
1211               endif
1213               CALL load_grid_info(TmpVarName, StartDate, vert_unit, level1(z),&
1214                    level2(z), fcst_secs, accum_period, wg_grid_id, &
1215                    projection, xsize, ysize, region_center_lat, &
1216                    region_center_lon, dx, dy, proj_central_lon, &
1217                    proj_center_flag, truelat1, truelat2, grib_tables, &
1218                    grid_info)
1219            
1220               CALL write_grib(grid_info, FileFd(DataHandle), data)
1221               CALL free_grid_info(grid_info)
1222            endif
1224         enddo
1225      endif
1226   endif
1228   deallocate(data, STAT = istat)
1229   deallocate(mold, STAT = istat)
1230   deallocate(tmpdata, STAT = istat)
1232   Status = WRF_NO_ERR
1234   call wrf_debug ( DEBUG , 'Leaving ext_gr1_write_field')
1236   RETURN
1237 END SUBROUTINE ext_gr1_write_field
1239 !*****************************************************************************
1241 SUBROUTINE ext_gr1_read_field ( DataHandle , DateStr , VarName , Field , &
1242      FieldType , Comm , IOComm, DomainDesc , MemoryOrder , Stagger ,     &
1243      DimNames , DomainStart , DomainEnd , MemoryStart , MemoryEnd ,      &
1244      PatchStart , PatchEnd ,  Status )
1246   USE gr1_data_info
1247   IMPLICIT NONE  
1248 #include "wrf_status_codes.h"
1249 #include "wrf_io_flags.h"
1250   INTEGER ,       INTENT(IN)    :: DataHandle 
1251   CHARACTER*(*) :: DateStr
1252   CHARACTER*(*) :: VarName
1253   CHARACTER (len=400) :: msg
1254   integer                       ,intent(inout)    :: FieldType
1255   integer                       ,intent(inout)    :: Comm
1256   integer                       ,intent(inout)    :: IOComm
1257   integer                       ,intent(inout)    :: DomainDesc
1258   character*(*)                 ,intent(inout)    :: MemoryOrder
1259   character*(*)                 ,intent(inout)    :: Stagger
1260   character*(*) , dimension (*) ,intent(inout)    :: DimNames
1261   integer ,dimension(*)         ,intent(inout)    :: DomainStart, DomainEnd
1262   integer ,dimension(*)         ,intent(inout)    :: MemoryStart, MemoryEnd
1263   integer ,dimension(*)         ,intent(inout)    :: PatchStart,  PatchEnd
1264   integer                       ,intent(out)      :: Status
1265   INTEGER                       ,intent(out)      :: Field(*)
1266   integer   :: ndim,x_start,x_end,y_start,y_end,z_start,z_end
1267   integer   :: zidx
1268   REAL, DIMENSION(:,:), POINTER :: data
1269   logical                     :: vert_stag
1270   logical                     :: soil_layers
1271   integer                     :: level1,level2
1273   integer                     :: parmid
1274   integer                     :: vert_unit
1275   integer                     :: grb_index
1276   integer                     :: numcols, numrows
1277   integer                     :: data_allocated
1278   integer                     :: istat
1279   integer                     :: tablenum
1280   integer                     :: di
1281   integer                     :: last_grb_index
1282   logical                     :: fraction
1284   call wrf_debug ( DEBUG , 'Entering ext_gr1_read_field')
1286   !
1287   ! Get dimensions of data.  
1288   ! Assume that the domain size in the input data is the same as the Domain 
1289   !     Size from the input arguments.
1290   !
1291   
1292   CALL get_dims(MemoryOrder,DomainStart,DomainEnd,ndim,x_start,x_end,y_start, &
1293        y_end,z_start,z_end) 
1295   !
1296   ! Get grib parameter id
1297   !
1298   CALL GET_GRIB_PARAM(grib_tables, VarName, center, subcenter, parmtbl, &
1299        tablenum, parmid)
1301   !
1302   ! Setup the vertical unit and levels
1303   !
1304   CALL get_vert_stag(VarName,Stagger,vert_stag)
1305   CALL get_soil_layers(VarName,soil_layers)
1307   !
1308   ! Loop over levels, grabbing data from each level, then assembling into a 
1309   !   3D array.
1310   !
1311   data_allocated = 0
1312   last_grb_index = -1
1313   do zidx = z_start,z_end
1314      
1315      IF ((VarName .eq. 'LANDUSEF') .or. (VarName .eq. 'SOILCBOT') .or. &
1316           (VarName .eq. 'SOILCTOP')) then
1317         fraction = .true.
1318      ELSE
1319         fraction = .false.
1320      END IF
1321      CALL gr1_get_levels(VarName,zidx,z_end-z_start,soil_layers,vert_stag,fraction, &
1322           .false., .false., vert_unit,level1,level2)
1323      
1324      CALL GET_GRIB_INDEX_VALIDTIME_GUESS(fileinfo(DataHandle)%fileindex(:), center, &
1325           subcenter, parmtbl, parmid,DateStr,vert_unit,level1, &
1326           level2, last_grb_index + 1, grb_index)
1327      if (grb_index < 0) then
1328         write(msg,*)'Field not found: parmid: ',VarName,parmid,DateStr, &
1329              vert_unit,level1,level2
1330         call wrf_debug (DEBUG , msg)
1331         cycle
1332      endif
1334      if (data_allocated .eq. 0) then
1335         CALL GET_SIZEOF_GRID(fileinfo(DataHandle)%fileindex(:),grb_index,numcols,numrows)
1336         allocate(data(z_start:z_end,1:numcols*numrows),stat=istat)
1337         data_allocated = 1
1338      endif
1340      CALL READ_GRIB(fileinfo(DataHandle)%fileindex(:), FileFd(DataHandle), grb_index, &
1341           data(zidx,:))
1343      !
1344      ! Transpose data into the order specified by MemoryOrder, setting only 
1345      !   entries within the memory dimensions
1346      !
1347      CALL get_dims(MemoryOrder, MemoryStart, MemoryEnd, ndim, x_start, x_end, &
1348           y_start, y_end,z_start,z_end)
1350      if(FieldType == WRF_DOUBLE) then
1351         di = 2
1352      else 
1353         di = 1
1354      endif
1356      ! 
1357      ! Here, we do any necessary conversions to the data.
1358      !
1359      ! The WRF executable (wrf.exe) expects perturbation potential
1360      !   temperature.  However, real.exe expects full potential T.
1361      ! So, if the program is WRF, subtract 300 from Potential Temperature 
1362      !   to get perturbation potential temperature.
1363      !
1364      if (VarName == 'T') then
1365         if ( &
1366              (InputProgramName .eq. 'REAL_EM') .or. &
1367              (InputProgramName .eq. 'IDEAL') .or. &
1368              (InputProgramName .eq. 'NDOWN_EM')) then
1369            data(zidx,:) = data(zidx,:) - 300
1370         endif
1371      endif
1373      CALL Transpose_grib(MemoryOrder, di, FieldType, Field, &
1374           MemoryStart(1), MemoryEnd(1), MemoryStart(2), MemoryEnd(2), &
1375           MemoryStart(3), MemoryEnd(3), &
1376           data(zidx,:), zidx, numrows, numcols)
1378      if (zidx .eq. z_end) then
1379         data_allocated = 0
1380         deallocate(data)
1381      endif
1383      last_grb_index = grb_index
1385   enddo
1387   Status = WRF_NO_ERR
1388   if (grb_index < 0) Status = WRF_WARN_VAR_NF
1389   call wrf_debug ( DEBUG , 'Leaving ext_gr1_read_field')
1391   RETURN
1392 END SUBROUTINE ext_gr1_read_field
1394 !*****************************************************************************
1396 SUBROUTINE ext_gr1_get_next_var ( DataHandle, VarName, Status )
1398   USE gr1_data_info
1399   IMPLICIT NONE
1400 #include "wrf_status_codes.h"
1401   INTEGER ,       INTENT(IN)  :: DataHandle
1402   CHARACTER*(*) :: VarName
1403   INTEGER ,       INTENT(OUT) :: Status
1405   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_next_var')
1407   call wrf_message ( 'WARNING: ext_gr1_get_next_var is not supported.')
1409   Status = WRF_WARN_NOOP
1411   RETURN
1412 END SUBROUTINE ext_gr1_get_next_var
1414 !*****************************************************************************
1416 subroutine ext_gr1_end_of_frame(DataHandle, Status)
1418   USE gr1_data_info
1419   implicit none
1420 #include "wrf_status_codes.h"
1421   integer               ,intent(in)     :: DataHandle
1422   integer               ,intent(out)    :: Status
1424   call wrf_debug ( DEBUG , 'Entering ext_gr1_end_of_frame')
1426   Status = WRF_WARN_NOOP
1428   return
1429 end subroutine ext_gr1_end_of_frame
1431 !*****************************************************************************
1433 SUBROUTINE ext_gr1_iosync ( DataHandle, Status )
1435   USE gr1_data_info  
1436   IMPLICIT NONE
1437 #include "wrf_status_codes.h"
1438   INTEGER ,       INTENT(IN)  :: DataHandle
1439   INTEGER ,       INTENT(OUT) :: Status
1441   call wrf_debug ( DEBUG , 'Entering ext_gr1_iosync')
1443   Status = WRF_NO_ERR
1444   if (DataHandle .GT. 0) then
1445      CALL flush_file(FileFd(DataHandle))
1446   else
1447      Status = WRF_WARN_TOO_MANY_FILES
1448   endif
1450   RETURN
1451 END SUBROUTINE ext_gr1_iosync
1453 !*****************************************************************************
1455 SUBROUTINE ext_gr1_inquire_filename ( DataHandle, FileName , FileStat, &
1456      Status )
1458   USE gr1_data_info
1459   IMPLICIT NONE
1460 #include "wrf_status_codes.h"
1461 #include "wrf_io_flags.h"
1462   INTEGER ,       INTENT(IN)  :: DataHandle
1463   CHARACTER*(*) :: FileName
1464   INTEGER ,       INTENT(OUT) :: FileStat
1465   INTEGER ,       INTENT(OUT) :: Status
1466   CHARACTER *80   SysDepInfo
1468   call wrf_debug ( DEBUG , 'Entering ext_gr1_inquire_filename')
1470   FileName = DataFile(DataHandle) 
1472   if ((DataHandle .ge. firstFileHandle) .and. &
1473        (DataHandle .le. maxFileHandles)) then
1474      FileStat = FileStatus(DataHandle)
1475   else
1476      FileStat = WRF_FILE_NOT_OPENED
1477   endif
1478   
1479   Status = WRF_NO_ERR
1481   RETURN
1482 END SUBROUTINE ext_gr1_inquire_filename
1484 !*****************************************************************************
1486 SUBROUTINE ext_gr1_get_var_info ( DataHandle , VarName , NDim , &
1487      MemoryOrder , Stagger , DomainStart , DomainEnd , Status )
1489   USE gr1_data_info
1490   IMPLICIT NONE
1491 #include "wrf_status_codes.h"
1492   integer               ,intent(in)     :: DataHandle
1493   character*(*)         ,intent(in)     :: VarName
1494   integer               ,intent(out)    :: NDim
1495   character*(*)         ,intent(out)    :: MemoryOrder
1496   character*(*)         ,intent(out)    :: Stagger
1497   integer ,dimension(*) ,intent(out)    :: DomainStart, DomainEnd
1498   integer               ,intent(out)    :: Status
1500   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_info')
1502   CALL wrf_message('ext_gr1_get_var_info not supported for grib version1 data')
1503   Status = WRF_NO_ERR
1505   RETURN
1506 END SUBROUTINE ext_gr1_get_var_info
1508 !*****************************************************************************
1510 SUBROUTINE ext_gr1_set_time ( DataHandle, DateStr, Status )
1512   USE gr1_data_info
1513   IMPLICIT NONE
1514 #include "wrf_status_codes.h"
1515   INTEGER ,       INTENT(IN)  :: DataHandle
1516   CHARACTER*(*) :: DateStr
1517   INTEGER ,       INTENT(OUT) :: Status
1518   integer       :: found_time
1519   integer       :: idx
1521   call wrf_debug ( DEBUG , 'Entering ext_gr1_set_time')
1523   found_time = 0
1524   do idx = 1,fileinfo(DataHandle)%NumberTimes
1525      if (fileinfo(DataHandle)%Times(idx) == DateStr) then
1526         found_time = 1
1527         fileinfo(DataHandle)%CurrentTime = idx
1528      endif
1529   enddo
1530   if (found_time == 0) then 
1531      Status = WRF_WARN_TIME_NF
1532   else
1533      Status = WRF_NO_ERR
1534   endif
1536   RETURN
1537 END SUBROUTINE ext_gr1_set_time
1539 !*****************************************************************************
1541 SUBROUTINE ext_gr1_get_next_time ( DataHandle, DateStr, Status )
1543   USE gr1_data_info
1544   IMPLICIT NONE
1545 #include "wrf_status_codes.h"
1546   INTEGER ,       INTENT(IN)  :: DataHandle
1547   CHARACTER*(*) , INTENT(OUT) :: DateStr
1548   INTEGER ,       INTENT(OUT) :: Status
1550   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_next_time')
1552   if (fileinfo(DataHandle)%CurrentTime == fileinfo(DataHandle)%NumberTimes) then
1553      Status = WRF_WARN_TIME_EOF
1554   else
1555      fileinfo(DataHandle)%CurrentTime = fileinfo(DataHandle)%CurrentTime + 1
1556      DateStr = fileinfo(DataHandle)%Times(fileinfo(DataHandle)%CurrentTime)
1557      Status = WRF_NO_ERR
1558   endif
1560   RETURN
1561 END SUBROUTINE ext_gr1_get_next_time
1563 !*****************************************************************************
1565 SUBROUTINE ext_gr1_get_previous_time ( DataHandle, DateStr, Status )
1567   USE gr1_data_info
1568   IMPLICIT NONE
1569 #include "wrf_status_codes.h"
1570   INTEGER ,       INTENT(IN)  :: DataHandle
1571   CHARACTER*(*) :: DateStr
1572   INTEGER ,       INTENT(OUT) :: Status
1574   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_previous_time')
1576   if (fileinfo(DataHandle)%CurrentTime <= 0) then
1577      Status = WRF_WARN_TIME_EOF
1578   else
1579      fileinfo(DataHandle)%CurrentTime = fileinfo(DataHandle)%CurrentTime - 1
1580      DateStr = fileinfo(DataHandle)%Times(fileinfo(DataHandle)%CurrentTime)
1581      Status = WRF_NO_ERR
1582   endif
1584   RETURN
1585 END SUBROUTINE ext_gr1_get_previous_time
1587 !******************************************************************************
1588 !* Start of get_var_ti_* routines
1589 !******************************************************************************
1591 SUBROUTINE ext_gr1_get_var_ti_real ( DataHandle,Element,  Varname, Data, &
1592      Count, Outcount, Status )
1594   USE gr1_data_info
1595   IMPLICIT NONE
1596 #include "wrf_status_codes.h"
1597   INTEGER ,       INTENT(IN)    :: DataHandle
1598   CHARACTER*(*) :: Element
1599   CHARACTER*(*) :: VarName 
1600   real ,          INTENT(OUT)   :: Data(*)
1601   INTEGER ,       INTENT(IN)    :: Count
1602   INTEGER ,       INTENT(OUT)   :: OutCount
1603   INTEGER ,       INTENT(OUT)   :: Status
1604   INTEGER          :: idx
1605   INTEGER          :: stat
1606   CHARACTER*(1000) :: VALUE
1608   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_ti_real')
1610   Status = WRF_NO_ERR
1611   
1612   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), "none", &
1613        Varname, Value, stat)
1614   if (stat /= 0) then
1615      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
1616      Status = WRF_WARN_VAR_NF
1617      RETURN
1618   endif
1620   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
1621   if (stat .ne. 0) then
1622      CALL wrf_message("Reading data from"//Value//"failed")
1623      Status = WRF_WARN_COUNT_TOO_LONG
1624      RETURN
1625   endif
1626   Outcount = idx
1628   RETURN
1629 END SUBROUTINE ext_gr1_get_var_ti_real 
1631 !*****************************************************************************
1633 SUBROUTINE ext_gr1_get_var_ti_real8 ( DataHandle,Element,  Varname, Data, &
1634      Count, Outcount, Status )
1636   USE gr1_data_info
1637   IMPLICIT NONE
1638 #include "wrf_status_codes.h"
1639   INTEGER ,       INTENT(IN)      :: DataHandle
1640   CHARACTER*(*) :: Element
1641   CHARACTER*(*) :: VarName 
1642   real*8 ,        INTENT(OUT)     :: Data(*)
1643   INTEGER ,       INTENT(IN)      :: Count
1644   INTEGER ,       INTENT(OUT)     :: OutCount
1645   INTEGER ,       INTENT(OUT)     :: Status
1646   INTEGER          :: idx
1647   INTEGER          :: stat
1648   CHARACTER*(1000) :: VALUE
1650   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_ti_real8')
1652   Status = WRF_NO_ERR
1653   
1654   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:),TRIM(Element),&
1655        "none",Varname,Value,stat)
1656   if (stat /= 0) then
1657      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
1658      Status = WRF_WARN_VAR_NF
1659      RETURN
1660   endif
1662   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
1663   if (stat .ne. 0) then
1664      CALL wrf_message("Reading data from"//Value//"failed")
1665      Status = WRF_WARN_COUNT_TOO_LONG
1666      RETURN
1667   endif
1668   Outcount = idx
1670   RETURN
1671 END SUBROUTINE ext_gr1_get_var_ti_real8 
1673 !*****************************************************************************
1675 SUBROUTINE ext_gr1_get_var_ti_double ( DataHandle,Element,  Varname, Data, &
1676      Count, Outcount, Status )
1677   USE gr1_data_info
1678   IMPLICIT NONE
1679 #include "wrf_status_codes.h"
1680   INTEGER ,       INTENT(IN)  :: DataHandle
1681   CHARACTER*(*) , INTENT(IN)  :: Element
1682   CHARACTER*(*) , INTENT(IN)  :: VarName
1683   real*8 ,            INTENT(OUT) :: Data(*)
1684   INTEGER ,       INTENT(IN)  :: Count
1685   INTEGER ,       INTENT(OUT)  :: OutCount
1686   INTEGER ,       INTENT(OUT) :: Status
1687   INTEGER          :: idx
1688   INTEGER          :: stat
1689   CHARACTER*(1000) :: VALUE
1691   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_ti_double')
1693   Status = WRF_NO_ERR
1694   
1695   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), &
1696        "none", Varname, &
1697        Value,stat)
1698   if (stat /= 0) then
1699      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
1700      Status = WRF_WARN_VAR_NF
1701      RETURN
1702   endif
1704   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
1705   if (stat .ne. 0) then
1706      CALL wrf_message("Reading data from"//Value//"failed")
1707      Status = WRF_WARN_COUNT_TOO_LONG
1708      RETURN
1709   endif
1710   Outcount = idx
1712   RETURN
1713 END SUBROUTINE ext_gr1_get_var_ti_double
1715 !*****************************************************************************
1717 SUBROUTINE ext_gr1_get_var_ti_integer ( DataHandle,Element,  Varname, Data, &
1718      Count, Outcount, Status )
1720   USE gr1_data_info
1721   IMPLICIT NONE
1722 #include "wrf_status_codes.h"
1723   INTEGER ,       INTENT(IN)       :: DataHandle
1724   CHARACTER*(*) :: Element
1725   CHARACTER*(*) :: VarName 
1726   integer ,       INTENT(OUT)      :: Data(*)
1727   INTEGER ,       INTENT(IN)       :: Count
1728   INTEGER ,       INTENT(OUT)      :: OutCount
1729   INTEGER ,       INTENT(OUT)      :: Status
1730   INTEGER          :: idx
1731   INTEGER          :: stat
1732   CHARACTER*(1000) :: VALUE
1734   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_ti_integer')
1736   Status = WRF_NO_ERR
1737   
1738   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), &
1739        "none", Varname, Value, stat)
1740   if (stat /= 0) then
1741      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
1742      Status = WRF_WARN_VAR_NF
1743      RETURN
1744   endif
1746   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
1747   if (stat .ne. 0) then
1748      CALL wrf_message("Reading data from"//Value//"failed")
1749      Status = WRF_WARN_COUNT_TOO_LONG
1750      RETURN
1751   endif
1752   Outcount = idx
1754   RETURN
1755 END SUBROUTINE ext_gr1_get_var_ti_integer 
1757 !*****************************************************************************
1759 SUBROUTINE ext_gr1_get_var_ti_logical ( DataHandle,Element,  Varname, Data, &
1760      Count, Outcount, Status )
1762   USE gr1_data_info
1763   IMPLICIT NONE
1764 #include "wrf_status_codes.h"
1765   INTEGER ,       INTENT(IN)       :: DataHandle
1766   CHARACTER*(*) :: Element
1767   CHARACTER*(*) :: VarName 
1768   logical ,       INTENT(OUT)      :: Data(*)
1769   INTEGER ,       INTENT(IN)       :: Count
1770   INTEGER ,       INTENT(OUT)      :: OutCount
1771   INTEGER ,       INTENT(OUT)      :: Status
1772   INTEGER          :: idx
1773   INTEGER          :: stat
1774   CHARACTER*(1000) :: VALUE
1776   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_ti_logical')
1778   Status = WRF_NO_ERR
1779   
1780   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), &
1781        "none", Varname, Value,stat)
1782   if (stat /= 0) then
1783      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
1784      Status = WRF_WARN_VAR_NF
1785      RETURN
1786   endif
1788   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
1789   if (stat .ne. 0) then
1790      CALL wrf_message("Reading data from"//Value//"failed")
1791      Status = WRF_WARN_COUNT_TOO_LONG
1792      RETURN
1793   endif
1794   Outcount = idx
1796   RETURN
1797 END SUBROUTINE ext_gr1_get_var_ti_logical 
1799 !*****************************************************************************
1801 SUBROUTINE ext_gr1_get_var_ti_char ( DataHandle,Element,  Varname, Data,  &
1802      Status )
1804   USE gr1_data_info
1805   IMPLICIT NONE
1806 #include "wrf_status_codes.h"
1807   INTEGER ,       INTENT(IN)  :: DataHandle
1808   CHARACTER*(*) :: Element
1809   CHARACTER*(*) :: VarName 
1810   CHARACTER*(*) :: Data
1811   INTEGER ,       INTENT(OUT) :: Status
1812   INTEGER       :: stat
1814   Status = WRF_NO_ERR
1815   
1816   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_ti_char')
1818   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), &
1819        "none", Varname, Data,stat)
1820   if (stat /= 0) then
1821      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
1822      Status = WRF_WARN_VAR_NF
1823      RETURN
1824   endif
1826   RETURN
1827 END SUBROUTINE ext_gr1_get_var_ti_char 
1829 !******************************************************************************
1830 !* End of get_var_ti_* routines
1831 !******************************************************************************
1834 !******************************************************************************
1835 !* Start of put_var_ti_* routines
1836 !******************************************************************************
1838 SUBROUTINE ext_gr1_put_var_ti_real ( DataHandle,Element,  Varname, Data, &
1839      Count,  Status )
1841   USE gr1_data_info
1842   IMPLICIT NONE
1843 #include "wrf_status_codes.h"
1844   INTEGER ,       INTENT(IN)  :: DataHandle
1845   CHARACTER*(*) :: Element
1846   CHARACTER*(*) :: VarName 
1847   real ,          INTENT(IN)  :: Data(*)
1848   INTEGER ,       INTENT(IN)  :: Count
1849   INTEGER ,       INTENT(OUT) :: Status
1850   CHARACTER(len=1000) :: tmpstr(1000)
1851   INTEGER             :: idx
1853   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_ti_real')
1855   if (committed(DataHandle)) then
1857      do idx = 1,Count
1858         write(tmpstr(idx),'(G17.10)')Data(idx)
1859      enddo
1861      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, Count, Status)
1863   endif
1865   RETURN
1866 END SUBROUTINE ext_gr1_put_var_ti_real 
1868 !*****************************************************************************
1870 SUBROUTINE ext_gr1_put_var_ti_double ( DataHandle,Element,  Varname, Data, &
1871      Count,  Status )
1872   USE gr1_data_info
1873   IMPLICIT NONE
1874 #include "wrf_status_codes.h"
1875   INTEGER ,       INTENT(IN)  :: DataHandle
1876   CHARACTER*(*) , INTENT(IN)  :: Element
1877   CHARACTER*(*) , INTENT(IN)  :: VarName
1878   real*8 ,            INTENT(IN) :: Data(*)
1879   INTEGER ,       INTENT(IN)  :: Count
1880   INTEGER ,       INTENT(OUT) :: Status
1881   CHARACTER(len=1000) :: tmpstr(1000)
1882   INTEGER             :: idx
1884   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_ti_double')
1886   if (committed(DataHandle)) then
1888      do idx = 1,Count
1889         write(tmpstr(idx),'(G17.10)')Data(idx)
1890      enddo
1891      
1892      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, Count, Status)
1893   endif
1895   RETURN
1896 END SUBROUTINE ext_gr1_put_var_ti_double
1898 !*****************************************************************************
1900 SUBROUTINE ext_gr1_put_var_ti_real8 ( DataHandle,Element,  Varname, Data, &
1901      Count,  Status )
1903   USE gr1_data_info
1904   IMPLICIT NONE
1905 #include "wrf_status_codes.h"
1906   INTEGER ,       INTENT(IN)  :: DataHandle
1907   CHARACTER*(*) :: Element
1908   CHARACTER*(*) :: VarName 
1909   real*8 ,        INTENT(IN)  :: Data(*)
1910   INTEGER ,       INTENT(IN)  :: Count
1911   INTEGER ,       INTENT(OUT) :: Status
1912   CHARACTER(len=1000) :: tmpstr(1000)
1913   INTEGER             :: idx
1915   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_ti_real8')
1917   if (committed(DataHandle)) then
1919      do idx = 1,Count
1920         write(tmpstr(idx),'(G17.10)')Data(idx)
1921      enddo
1922      
1923      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, Count, Status)
1924   endif
1926   RETURN
1927 END SUBROUTINE ext_gr1_put_var_ti_real8 
1929 !*****************************************************************************
1931 SUBROUTINE ext_gr1_put_var_ti_integer ( DataHandle,Element,  Varname, Data, &
1932      Count,  Status )
1934   USE gr1_data_info
1935   IMPLICIT NONE
1936 #include "wrf_status_codes.h"
1937   INTEGER ,       INTENT(IN)  :: DataHandle
1938   CHARACTER*(*) :: Element
1939   CHARACTER*(*) :: VarName 
1940   integer ,       INTENT(IN)  :: Data(*)
1941   INTEGER ,       INTENT(IN)  :: Count
1942   INTEGER ,       INTENT(OUT) :: Status
1943   CHARACTER(len=1000) :: tmpstr(1000)
1944   INTEGER             :: idx
1946   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_ti_integer')
1948   if (committed(DataHandle)) then
1950      do idx = 1,Count
1951         write(tmpstr(idx),'(G17.10)')Data(idx)
1952      enddo
1953      
1954      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, Count, Status)
1955   endif
1957   RETURN
1958 END SUBROUTINE ext_gr1_put_var_ti_integer 
1960 !*****************************************************************************
1962 SUBROUTINE ext_gr1_put_var_ti_logical ( DataHandle,Element,  Varname, Data, &
1963      Count,  Status )
1965   USE gr1_data_info
1966   IMPLICIT NONE
1967 #include "wrf_status_codes.h"
1968   INTEGER ,       INTENT(IN)  :: DataHandle
1969   CHARACTER*(*) :: Element
1970   CHARACTER*(*) :: VarName 
1971   logical ,       INTENT(IN)  :: Data(*)
1972   INTEGER ,       INTENT(IN)  :: Count
1973   INTEGER ,       INTENT(OUT) :: Status
1974   CHARACTER(len=1000) :: tmpstr(1000)
1975   INTEGER             :: idx
1977   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_ti_logical')
1979   if (committed(DataHandle)) then
1981      do idx = 1,Count
1982         write(tmpstr(idx),'(G17.10)')Data(idx)
1983      enddo
1984      
1985      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, Count, Status)
1987   endif
1989 RETURN
1990 END SUBROUTINE ext_gr1_put_var_ti_logical 
1992 !*****************************************************************************
1994 SUBROUTINE ext_gr1_put_var_ti_char ( DataHandle,Element,  Varname, Data,  &
1995      Status )
1997   USE gr1_data_info
1998   IMPLICIT NONE
1999 #include "wrf_status_codes.h"
2000   INTEGER ,       INTENT(IN)  :: DataHandle
2001   CHARACTER(len=*) :: Element
2002   CHARACTER(len=*) :: VarName 
2003   CHARACTER(len=*) :: Data
2004   INTEGER ,       INTENT(OUT) :: Status
2005   REAL dummy
2006   INTEGER                     :: Count
2007   CHARACTER(len=1000) :: tmpstr(1)
2008   INTEGER             :: idx
2010   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_ti_char')
2012   if (committed(DataHandle)) then
2014      write(tmpstr(1),*)trim(Data)
2016      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, 1, Status)
2018   endif
2020   RETURN
2021 END SUBROUTINE ext_gr1_put_var_ti_char 
2023 !******************************************************************************
2024 !* End of put_var_ti_* routines
2025 !******************************************************************************
2027 !******************************************************************************
2028 !* Start of get_var_td_* routines
2029 !******************************************************************************
2031 SUBROUTINE ext_gr1_get_var_td_double ( DataHandle,Element,  DateStr, &
2032      Varname, Data, Count, Outcount, Status )
2033   USE gr1_data_info
2034   IMPLICIT NONE
2035 #include "wrf_status_codes.h"
2036   INTEGER ,       INTENT(IN)  :: DataHandle
2037   CHARACTER*(*) , INTENT(IN)  :: Element
2038   CHARACTER*(*) , INTENT(IN)  :: DateStr
2039   CHARACTER*(*) , INTENT(IN)  :: VarName
2040   real*8 ,            INTENT(OUT) :: Data(*)
2041   INTEGER ,       INTENT(IN)  :: Count
2042   INTEGER ,       INTENT(OUT)  :: OutCount
2043   INTEGER ,       INTENT(OUT) :: Status
2044   INTEGER          :: idx
2045   INTEGER          :: stat
2046   CHARACTER*(1000) :: VALUE
2048   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_td_double')
2050   Status = WRF_NO_ERR
2051   
2052   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:),TRIM(Element),DateStr,&
2053        Varname,Value,stat)
2054   if (stat /= 0) then
2055      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2056      Status = WRF_WARN_VAR_NF
2057      RETURN
2058   endif
2060   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2061   if (stat .ne. 0) then
2062      CALL wrf_message("Reading data from"//Value//"failed")
2063      Status = WRF_WARN_COUNT_TOO_LONG
2064      RETURN
2065   endif
2066   Outcount = idx
2068 RETURN
2069 END SUBROUTINE ext_gr1_get_var_td_double
2071 !*****************************************************************************
2073 SUBROUTINE ext_gr1_get_var_td_real ( DataHandle,Element,  DateStr,Varname, &
2074      Data, Count, Outcount, Status )
2076   USE gr1_data_info
2077   IMPLICIT NONE
2078 #include "wrf_status_codes.h"
2079   INTEGER ,       INTENT(IN)  :: DataHandle
2080   CHARACTER*(*) :: Element
2081   CHARACTER*(*) :: DateStr
2082   CHARACTER*(*) :: VarName 
2083   real ,          INTENT(OUT) :: Data(*)
2084   INTEGER ,       INTENT(IN)  :: Count
2085   INTEGER ,       INTENT(OUT) :: OutCount
2086   INTEGER ,       INTENT(OUT) :: Status
2087   INTEGER          :: idx
2088   INTEGER          :: stat
2089   CHARACTER*(1000) :: VALUE
2091   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_td_real')
2093   Status = WRF_NO_ERR
2094   
2095   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), DateStr, &
2096        Varname, Value, stat)
2097   if (stat /= 0) then
2098      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2099      Status = WRF_WARN_VAR_NF
2100      RETURN
2101   endif
2103   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2104   if (stat .ne. 0) then
2105      CALL wrf_message("Reading data from"//Value//"failed")
2106      Status = WRF_WARN_COUNT_TOO_LONG
2107      RETURN
2108   endif
2109   Outcount = idx
2111   RETURN
2112 END SUBROUTINE ext_gr1_get_var_td_real 
2114 !*****************************************************************************
2116 SUBROUTINE ext_gr1_get_var_td_real8 ( DataHandle,Element,  DateStr,Varname, &
2117      Data, Count, Outcount, Status )
2119   USE gr1_data_info
2120   IMPLICIT NONE
2121 #include "wrf_status_codes.h"
2122   INTEGER ,       INTENT(IN)  :: DataHandle
2123   CHARACTER*(*) :: Element
2124   CHARACTER*(*) :: DateStr
2125   CHARACTER*(*) :: VarName 
2126   real*8 ,        INTENT(OUT) :: Data(*)
2127   INTEGER ,       INTENT(IN)  :: Count
2128   INTEGER ,       INTENT(OUT) :: OutCount
2129   INTEGER ,       INTENT(OUT) :: Status
2130   INTEGER          :: idx
2131   INTEGER          :: stat
2132   CHARACTER*(1000) :: VALUE
2134   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_td_real8')
2136   Status = WRF_NO_ERR
2137   
2138   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:),TRIM(Element),DateStr,&
2139        Varname,Value,stat)
2140   if (stat /= 0) then
2141      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2142      Status = WRF_WARN_VAR_NF
2143      RETURN
2144   endif
2146   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2147   if (stat .ne. 0) then
2148      CALL wrf_message("Reading data from"//Value//"failed")
2149      Status = WRF_WARN_COUNT_TOO_LONG
2150      RETURN
2151   endif
2152   Outcount = idx
2154   RETURN
2155 END SUBROUTINE ext_gr1_get_var_td_real8 
2157 !*****************************************************************************
2159 SUBROUTINE ext_gr1_get_var_td_integer ( DataHandle,Element,  DateStr,Varname, &
2160      Data, Count, Outcount, Status )
2162   USE gr1_data_info
2163   IMPLICIT NONE
2164 #include "wrf_status_codes.h"
2165   INTEGER ,       INTENT(IN)  :: DataHandle
2166   CHARACTER*(*) :: Element
2167   CHARACTER*(*) :: DateStr
2168   CHARACTER*(*) :: VarName 
2169   integer ,       INTENT(OUT) :: Data(*)
2170   INTEGER ,       INTENT(IN)  :: Count
2171   INTEGER ,       INTENT(OUT) :: OutCount
2172   INTEGER ,       INTENT(OUT) :: Status
2173   INTEGER          :: idx
2174   INTEGER          :: stat
2175   CHARACTER*(1000) :: VALUE
2177   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_td_integer')
2179   Status = WRF_NO_ERR
2180   
2181   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), DateStr, &
2182        Varname, Value,stat)
2183   if (stat /= 0) then
2184      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2185      Status = WRF_WARN_VAR_NF
2186      RETURN
2187   endif
2189   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2190   if (stat .ne. 0) then
2191      CALL wrf_message("Reading data from"//Value//"failed")
2192      Status = WRF_WARN_COUNT_TOO_LONG
2193      RETURN
2194   endif
2195   Outcount = idx
2197   RETURN
2198 END SUBROUTINE ext_gr1_get_var_td_integer 
2200 !*****************************************************************************
2202 SUBROUTINE ext_gr1_get_var_td_logical ( DataHandle,Element,  DateStr,Varname, &
2203      Data, Count, Outcount, Status )
2204   
2205   USE gr1_data_info
2206   IMPLICIT NONE
2207 #include "wrf_status_codes.h"
2208   INTEGER ,       INTENT(IN)  :: DataHandle
2209   CHARACTER*(*) :: Element
2210   CHARACTER*(*) :: DateStr
2211   CHARACTER*(*) :: VarName 
2212   logical ,       INTENT(OUT) :: Data(*)
2213   INTEGER ,       INTENT(IN)  :: Count
2214   INTEGER ,       INTENT(OUT) :: OutCount
2215   INTEGER ,       INTENT(OUT) :: Status
2216   INTEGER          :: idx
2217   INTEGER          :: stat
2218   CHARACTER*(1000) :: VALUE
2220   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_td_logical')
2222   Status = WRF_NO_ERR
2223   
2224   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), DateStr, &
2225        Varname, Value,stat)
2226   if (stat /= 0) then
2227      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2228      Status = WRF_WARN_VAR_NF
2229      RETURN
2230   endif
2232   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2233   if (stat .ne. 0) then
2234      CALL wrf_message("Reading data from"//Value//"failed")
2235      Status = WRF_WARN_COUNT_TOO_LONG
2236      RETURN
2237   endif
2238   Outcount = idx
2240   RETURN
2241 END SUBROUTINE ext_gr1_get_var_td_logical 
2243 !*****************************************************************************
2245 SUBROUTINE ext_gr1_get_var_td_char ( DataHandle,Element,  DateStr,Varname, &
2246      Data,  Status )
2248   USE gr1_data_info
2249   IMPLICIT NONE
2250 #include "wrf_status_codes.h"
2251   INTEGER ,       INTENT(IN)  :: DataHandle
2252   CHARACTER*(*) :: Element
2253   CHARACTER*(*) :: DateStr
2254   CHARACTER*(*) :: VarName 
2255   CHARACTER*(*) :: Data
2256   INTEGER ,       INTENT(OUT) :: Status
2257   INTEGER       :: stat
2259   Status = WRF_NO_ERR
2260   
2261   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_td_char')
2263   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), DateStr, &
2264        Varname, Data,stat)
2265   if (stat /= 0) then
2266      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2267      Status = WRF_WARN_VAR_NF
2268      RETURN
2269   endif
2271   RETURN
2272 END SUBROUTINE ext_gr1_get_var_td_char 
2274 !******************************************************************************
2275 !* End of get_var_td_* routines
2276 !******************************************************************************
2278 !******************************************************************************
2279 !* Start of put_var_td_* routines
2280 !******************************************************************************
2282 SUBROUTINE ext_gr1_put_var_td_double ( DataHandle, Element, DateStr, Varname, &
2283      Data, Count,  Status )
2284   USE gr1_data_info
2285   IMPLICIT NONE
2286 #include "wrf_status_codes.h"
2287   INTEGER ,       INTENT(IN)  :: DataHandle
2288   CHARACTER*(*) , INTENT(IN)  :: Element
2289   CHARACTER*(*) , INTENT(IN)  :: DateStr
2290   CHARACTER*(*) , INTENT(IN)  :: VarName
2291   real*8 ,            INTENT(IN) :: Data(*)
2292   INTEGER ,       INTENT(IN)  :: Count
2293   INTEGER ,       INTENT(OUT) :: Status
2294   CHARACTER(len=1000) :: tmpstr(1000)
2295   INTEGER             :: idx
2297   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_td_double')
2300   if (committed(DataHandle)) then
2302      do idx = 1,Count
2303         write(tmpstr(idx),'(G17.10)')Data(idx)
2304      enddo
2306      CALL gr1_build_string (td_output(DataHandle), &
2307           Varname//';'//DateStr//';'//Element, tmpstr, Count, Status)
2309   endif
2311 RETURN
2312 END SUBROUTINE ext_gr1_put_var_td_double
2314 !*****************************************************************************
2316 SUBROUTINE ext_gr1_put_var_td_integer ( DataHandle,Element,  DateStr, &
2317      Varname, Data, Count,  Status )
2319   USE gr1_data_info
2320   IMPLICIT NONE
2321 #include "wrf_status_codes.h"
2322   INTEGER ,       INTENT(IN)  :: DataHandle
2323   CHARACTER*(*) :: Element
2324   CHARACTER*(*) :: DateStr
2325   CHARACTER*(*) :: VarName 
2326   integer ,       INTENT(IN)  :: Data(*)
2327   INTEGER ,       INTENT(IN)  :: Count
2328   INTEGER ,       INTENT(OUT) :: Status
2329   CHARACTER(len=1000) :: tmpstr(1000)
2330   INTEGER             :: idx
2332   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_td_integer')
2334   if (committed(DataHandle)) then
2336      do idx = 1,Count
2337         write(tmpstr(idx),'(G17.10)')Data(idx)
2338      enddo
2339      
2340      CALL gr1_build_string (td_output(DataHandle), &
2341           Varname//';'//DateStr//';'//Element, tmpstr, Count, Status)
2343   endif
2345 RETURN
2346 END SUBROUTINE ext_gr1_put_var_td_integer 
2348 !*****************************************************************************
2350 SUBROUTINE ext_gr1_put_var_td_real ( DataHandle,Element,  DateStr,Varname, &
2351      Data, Count,  Status )
2353   USE gr1_data_info
2354   IMPLICIT NONE
2355 #include "wrf_status_codes.h"
2356   INTEGER ,       INTENT(IN)  :: DataHandle
2357   CHARACTER*(*) :: Element
2358   CHARACTER*(*) :: DateStr
2359   CHARACTER*(*) :: VarName 
2360   real ,          INTENT(IN)  :: Data(*)
2361   INTEGER ,       INTENT(IN)  :: Count
2362   INTEGER ,       INTENT(OUT) :: Status
2363   CHARACTER(len=1000) :: tmpstr(1000)
2364   INTEGER             :: idx
2366   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_td_real')
2368   if (committed(DataHandle)) then
2370      do idx = 1,Count
2371         write(tmpstr(idx),'(G17.10)')Data(idx)
2372      enddo
2373      
2374      CALL gr1_build_string (td_output(DataHandle), &
2375           Varname//';'//DateStr//';'//Element, tmpstr, Count, Status)
2377   endif
2379   RETURN
2380 END SUBROUTINE ext_gr1_put_var_td_real 
2382 !*****************************************************************************
2384 SUBROUTINE ext_gr1_put_var_td_real8 ( DataHandle,Element,  DateStr,Varname, &
2385      Data, Count,  Status )
2387   USE gr1_data_info
2388   IMPLICIT NONE
2389 #include "wrf_status_codes.h"
2390   INTEGER ,       INTENT(IN)  :: DataHandle
2391   CHARACTER*(*) :: Element
2392   CHARACTER*(*) :: DateStr
2393   CHARACTER*(*) :: VarName 
2394   real*8 ,        INTENT(IN)  :: Data(*)
2395   INTEGER ,       INTENT(IN)  :: Count
2396   INTEGER ,       INTENT(OUT) :: Status
2397   CHARACTER(len=1000) :: tmpstr(1000)
2398   INTEGER             :: idx
2400   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_td_real8')
2402   if (committed(DataHandle)) then
2403      do idx = 1,Count
2404         write(tmpstr(idx),'(G17.10)')Data(idx)
2405      enddo
2406      
2407      CALL gr1_build_string (td_output(DataHandle), &
2408           Varname//';'//DateStr//';'//Element, tmpstr, Count, Status)
2409   endif
2411   RETURN
2412 END SUBROUTINE ext_gr1_put_var_td_real8 
2414 !*****************************************************************************
2416 SUBROUTINE ext_gr1_put_var_td_logical ( DataHandle,Element,  DateStr, &
2417      Varname, Data, Count,  Status )
2419   USE gr1_data_info
2420   IMPLICIT NONE
2421 #include "wrf_status_codes.h"
2422   INTEGER ,       INTENT(IN)  :: DataHandle
2423   CHARACTER*(*) :: Element
2424   CHARACTER*(*) :: DateStr
2425   CHARACTER*(*) :: VarName 
2426   logical ,       INTENT(IN)  :: Data(*)
2427   INTEGER ,       INTENT(IN)  :: Count
2428   INTEGER ,       INTENT(OUT) :: Status
2429   CHARACTER(len=1000) :: tmpstr(1000)
2430   INTEGER             :: idx
2432   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_td_logical')
2434   if (committed(DataHandle)) then
2436      do idx = 1,Count
2437         write(tmpstr(idx),'(G17.10)')Data(idx)
2438      enddo
2440      CALL gr1_build_string (td_output(DataHandle), &
2441           Varname//';'//DateStr//';'//Element, tmpstr, Count, Status)
2443   endif
2445   RETURN
2446 END SUBROUTINE ext_gr1_put_var_td_logical 
2448 !*****************************************************************************
2450 SUBROUTINE ext_gr1_put_var_td_char ( DataHandle,Element,  DateStr,Varname, &
2451      Data,  Status )
2453   USE gr1_data_info
2454   IMPLICIT NONE
2455 #include "wrf_status_codes.h"
2456   INTEGER ,       INTENT(IN)  :: DataHandle
2457   CHARACTER*(*) :: Element
2458   CHARACTER*(*) :: DateStr
2459   CHARACTER*(*) :: VarName 
2460   CHARACTER*(*) :: Data
2461   INTEGER ,       INTENT(OUT) :: Status
2462   CHARACTER(len=1000) :: tmpstr
2463   INTEGER             :: idx
2465   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_td_char')
2467   if (committed(DataHandle)) then
2469     
2470      DO idx=1,LEN(Data)
2471         tmpstr(idx:idx)=Data(idx:idx)
2472      END DO
2473      DO idx=LEN(Data)+1,1000
2474         tmpstr(idx:idx)=' '
2475      END DO
2477      CALL gr1_build_string (td_output(DataHandle), &
2478           Varname//';'//DateStr//';'//Element, tmpstr, 1, Status)
2480   endif
2482   RETURN
2483 END SUBROUTINE ext_gr1_put_var_td_char 
2485 !******************************************************************************
2486 !* End of put_var_td_* routines
2487 !******************************************************************************
2490 !******************************************************************************
2491 !* Start of get_dom_ti_* routines
2492 !******************************************************************************
2494 SUBROUTINE ext_gr1_get_dom_ti_real ( DataHandle,Element,   Data, Count, &
2495      Outcount, Status )
2497   USE gr1_data_info
2498   IMPLICIT NONE
2499 #include "wrf_status_codes.h"
2500   INTEGER ,       INTENT(IN)  :: DataHandle
2501   CHARACTER*(*) :: Element
2502   real ,          INTENT(OUT) :: Data(*)
2503   INTEGER ,       INTENT(IN)  :: Count
2504   INTEGER ,       INTENT(OUT) :: Outcount
2505   INTEGER ,       INTENT(OUT) :: Status
2506   INTEGER          :: idx
2507   INTEGER          :: stat
2508   CHARACTER*(1000) :: VALUE
2510   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_ti_real')
2512   Status = WRF_NO_ERR
2513   
2514   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), "none", &
2515        "none", Value,stat)
2516   if (stat /= 0) then
2517      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2518      Status = WRF_WARN_VAR_NF
2519      RETURN
2520   endif
2522   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2523   if (stat .ne. 0) then
2524      CALL wrf_message("Reading data from"//Value//"failed")
2525      Status = WRF_WARN_COUNT_TOO_LONG
2526      RETURN
2527   endif
2528   Outcount = idx
2530   RETURN
2531 END SUBROUTINE ext_gr1_get_dom_ti_real 
2533 !*****************************************************************************
2535 SUBROUTINE ext_gr1_get_dom_ti_real8 ( DataHandle,Element,   Data, Count, &
2536      Outcount, Status )
2538   USE gr1_data_info
2539   IMPLICIT NONE
2540 #include "wrf_status_codes.h"
2541   INTEGER ,       INTENT(IN)  :: DataHandle
2542   CHARACTER*(*) :: Element
2543   real*8 ,        INTENT(OUT) :: Data(*)
2544   INTEGER ,       INTENT(IN)  :: Count
2545   INTEGER ,       INTENT(OUT) :: OutCount
2546   INTEGER ,       INTENT(OUT) :: Status
2547   INTEGER          :: idx
2548   INTEGER          :: stat
2549   CHARACTER*(1000) :: VALUE
2551   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_ti_real8')
2553   Status = WRF_NO_ERR
2554   
2555   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), "none", &
2556        "none", Value,stat)
2557   if (stat /= 0) then
2558      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2559      Status = WRF_WARN_VAR_NF
2560      RETURN
2561   endif
2563   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2564   if (stat .ne. 0) then
2565      CALL wrf_message("Reading data from"//Value//"failed")
2566      Status = WRF_WARN_COUNT_TOO_LONG
2567      RETURN
2568   endif
2569   Outcount = idx
2571   RETURN
2572 END SUBROUTINE ext_gr1_get_dom_ti_real8 
2574 !*****************************************************************************
2576 SUBROUTINE ext_gr1_get_dom_ti_integer ( DataHandle,Element,   Data, Count, &
2577      Outcount, Status )
2579   USE gr1_data_info
2580   IMPLICIT NONE
2581 #include "wrf_status_codes.h"
2582   INTEGER ,       INTENT(IN)  :: DataHandle
2583   CHARACTER*(*) :: Element
2584   integer ,       INTENT(OUT) :: Data(*)
2585   INTEGER ,       INTENT(IN)  :: Count
2586   INTEGER ,       INTENT(OUT) :: OutCount
2587   INTEGER ,       INTENT(OUT) :: Status
2588   INTEGER          :: idx
2589   INTEGER          :: stat
2590   CHARACTER*(1000) :: VALUE
2591   
2592   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_ti_integer Element: '//Element)
2594   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), "none", &
2595        "none", Value,stat)
2596   if (stat /= 0) then
2597      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2598      Status = WRF_WARN_VAR_NF
2599      RETURN
2600   endif
2602   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2603   if (stat .ne. 0) then
2604      CALL wrf_message("Reading data from"//Value//"failed")
2605      Status = WRF_WARN_COUNT_TOO_LONG
2606      RETURN
2607   endif
2608   Outcount = Count
2610   RETURN
2611 END SUBROUTINE ext_gr1_get_dom_ti_integer 
2613 !*****************************************************************************
2615 SUBROUTINE ext_gr1_get_dom_ti_logical ( DataHandle,Element,   Data, Count, &
2616      Outcount, Status )
2618   USE gr1_data_info
2619   IMPLICIT NONE
2620 #include "wrf_status_codes.h"
2621   INTEGER ,       INTENT(IN)  :: DataHandle
2622   CHARACTER*(*) :: Element
2623   logical ,       INTENT(OUT) :: Data(*)
2624   INTEGER ,       INTENT(IN)  :: Count
2625   INTEGER ,       INTENT(OUT) :: OutCount
2626   INTEGER ,       INTENT(OUT) :: Status
2627   INTEGER          :: idx
2628   INTEGER          :: stat
2629   CHARACTER*(1000) :: VALUE
2631   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_ti_logical')
2633   Status = WRF_NO_ERR
2634   
2635   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), "none", &
2636        "none", Value,stat)
2637   if (stat /= 0) then
2638      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2639      Status = WRF_WARN_VAR_NF
2640      RETURN
2641   endif
2643   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2644   if (stat .ne. 0) then
2645      CALL wrf_message("Reading data from"//Value//"failed")
2646      Status = WRF_WARN_COUNT_TOO_LONG
2647      RETURN
2648   endif
2649   Outcount = idx
2651   RETURN
2652 END SUBROUTINE ext_gr1_get_dom_ti_logical 
2654 !*****************************************************************************
2656 SUBROUTINE ext_gr1_get_dom_ti_char ( DataHandle,Element,   Data,  Status )
2658   USE gr1_data_info
2659   IMPLICIT NONE
2660 #include "wrf_status_codes.h"
2661   INTEGER ,       INTENT(IN)  :: DataHandle
2662   CHARACTER*(*) :: Element
2663   CHARACTER*(*) :: Data
2664   INTEGER ,       INTENT(OUT) :: Status
2665   INTEGER       :: stat
2666   INTEGER       :: endchar
2668   Status = WRF_NO_ERR
2669   
2670   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_ti_char')
2672   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), "none", &
2673        "none", Data, stat)
2674   if (stat /= 0) then
2675      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2676      Status = WRF_WARN_VAR_NF
2677      RETURN
2678   endif
2680   RETURN
2681 END SUBROUTINE ext_gr1_get_dom_ti_char 
2683 !*****************************************************************************
2685 SUBROUTINE ext_gr1_get_dom_ti_double ( DataHandle,Element,   Data, Count, &
2686      Outcount, Status )
2687   USE gr1_data_info
2688   IMPLICIT NONE
2689 #include "wrf_status_codes.h"
2690   INTEGER ,       INTENT(IN)  :: DataHandle
2691   CHARACTER*(*) , INTENT(IN)  :: Element
2692   real*8 ,            INTENT(OUT) :: Data(*)
2693   INTEGER ,       INTENT(IN)  :: Count
2694   INTEGER ,       INTENT(OUT)  :: OutCount
2695   INTEGER ,       INTENT(OUT) :: Status
2696   INTEGER          :: idx
2697   INTEGER          :: stat
2698   CHARACTER*(1000) :: VALUE
2700   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_ti_double')
2702   Status = WRF_NO_ERR
2703   
2704   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), "none", &
2705        "none", Value, stat)
2706   if (stat /= 0) then
2707      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2708      Status = WRF_WARN_VAR_NF
2709      RETURN
2710   endif
2712   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2713   if (stat .ne. 0) then
2714      CALL wrf_message("Reading data from"//Value//"failed")
2715      Status = WRF_WARN_COUNT_TOO_LONG
2716      RETURN
2717   endif
2718   Outcount = idx
2720 RETURN
2721 END SUBROUTINE ext_gr1_get_dom_ti_double
2723 !******************************************************************************
2724 !* End of get_dom_ti_* routines
2725 !******************************************************************************
2728 !******************************************************************************
2729 !* Start of put_dom_ti_* routines
2730 !******************************************************************************
2732 SUBROUTINE ext_gr1_put_dom_ti_real ( DataHandle,Element,   Data, Count,  &
2733      Status )
2735   USE gr1_data_info
2736   IMPLICIT NONE
2737 #include "wrf_status_codes.h"
2738   INTEGER ,       INTENT(IN)  :: DataHandle
2739   CHARACTER*(*) :: Element
2740   real ,          INTENT(IN)  :: Data(*)
2741   INTEGER ,       INTENT(IN)  :: Count
2742   INTEGER ,       INTENT(OUT) :: Status
2743   REAL dummy
2744   CHARACTER(len=1000) :: tmpstr(1000)
2745   character(len=2)    :: lf
2746   integer             :: idx
2748   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_ti_real')
2750   if (Element .eq. 'DX') then
2751      dx = Data(1)/1000.
2752   endif
2753   if (Element .eq. 'DY') then
2754      dy = Data(1)/1000.
2755   endif
2756   if (Element .eq. 'CEN_LAT') then
2757      center_lat = Data(1)
2758   endif
2759   if (Element .eq. 'CEN_LON') then
2760      center_lon = Data(1)
2761   endif  
2762   if (Element .eq. 'TRUELAT1') then
2763      truelat1 = Data(1)
2764   endif
2765   if (Element .eq. 'TRUELAT2') then
2766      truelat2 = Data(1)
2767   endif
2768   if (Element == 'STAND_LON') then
2769      proj_central_lon = Data(1)
2770   endif
2771   if (Element == 'DT') then
2772      timestep = Data(1)
2773   endif
2775   if (committed(DataHandle)) then
2777      do idx = 1,Count
2778         write(tmpstr(idx),'(G17.10)')Data(idx)
2779      enddo
2780      
2781      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, Count, Status)
2783   endif
2785   RETURN
2786 END SUBROUTINE ext_gr1_put_dom_ti_real 
2788 !*****************************************************************************
2790 SUBROUTINE ext_gr1_put_dom_ti_real8 ( DataHandle,Element,   Data, Count,  &
2791      Status )
2793   USE gr1_data_info
2794   IMPLICIT NONE
2795 #include "wrf_status_codes.h"
2796   INTEGER ,       INTENT(IN)  :: DataHandle
2797   CHARACTER*(*) :: Element
2798   real*8 ,        INTENT(IN)  :: Data(*)
2799   INTEGER ,       INTENT(IN)  :: Count
2800   INTEGER ,       INTENT(OUT) :: Status
2801   CHARACTER(len=1000) :: tmpstr(1000)
2802   INTEGER             :: idx
2804   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_ti_real8')
2806   if (committed(DataHandle)) then
2808      do idx = 1,Count
2809         write(tmpstr(idx),'(G17.10)')Data(idx)
2810      enddo
2811      
2812      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, Count, Status)
2814   endif
2816   RETURN
2817 END SUBROUTINE ext_gr1_put_dom_ti_real8 
2819 !*****************************************************************************
2821 SUBROUTINE ext_gr1_put_dom_ti_integer ( DataHandle,Element,   Data, Count,  &
2822      Status )
2824   USE gr1_data_info
2825   IMPLICIT NONE
2826 #include "wrf_status_codes.h"
2827   INTEGER ,       INTENT(IN)  :: DataHandle
2828   CHARACTER*(*) :: Element
2829   INTEGER ,       INTENT(IN)  :: Data(*)
2830   INTEGER ,       INTENT(IN)  :: Count
2831   INTEGER ,       INTENT(OUT) :: Status
2832   REAL dummy
2833   CHARACTER(len=1000) :: tmpstr(1000)
2834   INTEGER             :: idx
2837   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_ti_integer')
2839   if (Element == 'WEST-EAST_GRID_DIMENSION') then
2840      full_xsize = Data(1)
2841   else if (Element == 'SOUTH-NORTH_GRID_DIMENSION') then
2842      full_ysize = Data(1)
2843   else if (Element == 'MAP_PROJ') then
2844      projection = Data(1)
2845   else if (Element == 'WG_GRID_ID') then
2846      wg_grid_id = Data(1)
2847   else if (Element == 'GRID_ID') then
2848      this_domain = Data(1)
2849   endif
2851   if (committed(DataHandle)) then
2853      do idx = 1,Count
2854         write(tmpstr(idx),'(G17.10)')Data(idx)
2855      enddo
2856      
2857      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, Count, Status)
2859   endif
2861   call wrf_debug ( DEBUG , 'Leaving ext_gr1_put_dom_ti_integer')
2863   RETURN
2864 END SUBROUTINE ext_gr1_put_dom_ti_integer 
2866 !*****************************************************************************
2868 SUBROUTINE ext_gr1_put_dom_ti_logical ( DataHandle,Element,   Data, Count,  &
2869      Status )
2871   USE gr1_data_info
2872   IMPLICIT NONE
2873 #include "wrf_status_codes.h"
2874   INTEGER ,       INTENT(IN)  :: DataHandle
2875   CHARACTER*(*) :: Element
2876   logical ,       INTENT(IN)  :: Data(*)
2877   INTEGER ,       INTENT(IN)  :: Count
2878   INTEGER ,       INTENT(OUT) :: Status
2879   CHARACTER(len=1000) :: tmpstr(1000)
2880   INTEGER             :: idx
2882   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_ti_logical')
2884   if (committed(DataHandle)) then
2886      do idx = 1,Count
2887         write(tmpstr(idx),'(G17.10)')Data(idx)
2888      enddo
2889      
2890      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, Count, Status)
2892   endif
2894   RETURN
2895 END SUBROUTINE ext_gr1_put_dom_ti_logical 
2897 !*****************************************************************************
2899 SUBROUTINE ext_gr1_put_dom_ti_char ( DataHandle,Element,   Data,  &
2900      Status )
2902   USE gr1_data_info
2903   IMPLICIT NONE
2904 #include "wrf_status_codes.h"
2905   INTEGER ,       INTENT(IN)  :: DataHandle
2906   CHARACTER*(*) :: Element
2907   CHARACTER*(*),     INTENT(IN)  :: Data
2908   INTEGER ,       INTENT(OUT) :: Status
2909   REAL dummy
2910   CHARACTER(len=1000) :: tmpstr(1000)
2912   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_ti_char')
2914   if (Element .eq. 'START_DATE') then
2915      StartDate = Data
2916   endif
2918   if (committed(DataHandle)) then
2920      write(tmpstr(1),*)trim(Data)
2921      
2922      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, 1, Status)
2924   endif
2926   RETURN
2927 END SUBROUTINE ext_gr1_put_dom_ti_char
2929 !*****************************************************************************
2931 SUBROUTINE ext_gr1_put_dom_ti_double ( DataHandle,Element, Data, Count, &
2932      Status )
2933   USE gr1_data_info
2934   IMPLICIT NONE
2935 #include "wrf_status_codes.h"
2936   INTEGER ,       INTENT(IN)  :: DataHandle
2937   CHARACTER*(*) , INTENT(IN)  :: Element
2938   real*8 ,            INTENT(IN) :: Data(*)
2939   INTEGER ,       INTENT(IN)  :: Count
2940   INTEGER ,       INTENT(OUT) :: Status
2941   CHARACTER(len=1000) :: tmpstr(1000)
2942   INTEGER             :: idx
2944   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_ti_double')
2946   if (committed(DataHandle)) then
2948      do idx = 1,Count
2949         write(tmpstr(idx),'(G17.10)')Data(idx)
2950      enddo
2952      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, Count, Status)
2954   endif
2955   
2956   RETURN
2957 END SUBROUTINE ext_gr1_put_dom_ti_double
2959 !******************************************************************************
2960 !* End of put_dom_ti_* routines
2961 !******************************************************************************
2964 !******************************************************************************
2965 !* Start of get_dom_td_* routines
2966 !******************************************************************************
2968 SUBROUTINE ext_gr1_get_dom_td_real ( DataHandle,Element, DateStr,  Data, &
2969      Count, Outcount, Status )
2971   USE gr1_data_info
2972   IMPLICIT NONE
2973 #include "wrf_status_codes.h"
2974   INTEGER ,       INTENT(IN)  :: DataHandle
2975   CHARACTER*(*) :: Element
2976   CHARACTER*(*) :: DateStr
2977   real ,          INTENT(OUT) :: Data(*)
2978   INTEGER ,       INTENT(IN)  :: Count
2979   INTEGER ,       INTENT(OUT) :: OutCount
2980   INTEGER ,       INTENT(OUT) :: Status
2981   INTEGER          :: idx
2982   INTEGER          :: stat
2983   CHARACTER*(1000) :: VALUE
2985   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_td_real')
2987   Status = WRF_NO_ERR
2988   
2989   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), DateStr, &
2990        "none", Value, stat)
2991   if (stat /= 0) then
2992      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2993      Status = WRF_WARN_VAR_NF
2994      RETURN
2995   endif
2997   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2998   if (stat .ne. 0) then
2999      CALL wrf_message("Reading data from"//Value//"failed")
3000      Status = WRF_WARN_COUNT_TOO_LONG
3001      RETURN
3002   endif
3003   Outcount = idx
3005   RETURN
3006 END SUBROUTINE ext_gr1_get_dom_td_real 
3008 !*****************************************************************************
3010 SUBROUTINE ext_gr1_get_dom_td_real8 ( DataHandle,Element, DateStr,  Data, &
3011      Count, Outcount, Status )
3013   USE gr1_data_info
3014   IMPLICIT NONE
3015 #include "wrf_status_codes.h"
3016   INTEGER ,       INTENT(IN)  :: DataHandle
3017   CHARACTER*(*) :: Element
3018   CHARACTER*(*) :: DateStr
3019   real*8 ,        INTENT(OUT) :: Data(*)
3020   INTEGER ,       INTENT(IN)  :: Count
3021   INTEGER ,       INTENT(OUT) :: OutCount
3022   INTEGER ,       INTENT(OUT) :: Status
3023   INTEGER          :: idx
3024   INTEGER          :: stat
3025   CHARACTER*(1000) :: VALUE
3027   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_td_real8')
3029   Status = WRF_NO_ERR
3030   
3031   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), DateStr, &
3032        "none", Value, stat)
3033   if (stat /= 0) then
3034      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
3035      Status = WRF_WARN_VAR_NF
3036      RETURN
3037   endif
3039   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
3040   if (stat .ne. 0) then
3041      CALL wrf_message("Reading data from"//Value//"failed")
3042      Status = WRF_WARN_COUNT_TOO_LONG
3043      RETURN
3044   endif
3045   Outcount = idx
3047   RETURN
3048 END SUBROUTINE ext_gr1_get_dom_td_real8 
3050 !*****************************************************************************
3052 SUBROUTINE ext_gr1_get_dom_td_integer ( DataHandle,Element, DateStr,  Data, &
3053      Count, Outcount, Status )
3055   USE gr1_data_info
3056   IMPLICIT NONE
3057 #include "wrf_status_codes.h"
3058   INTEGER ,       INTENT(IN)  :: DataHandle
3059   CHARACTER*(*) :: Element
3060   CHARACTER*(*) :: DateStr
3061   integer ,       INTENT(OUT) :: Data(*)
3062   INTEGER ,       INTENT(IN)  :: Count
3063   INTEGER ,       INTENT(OUT) :: OutCount
3064   INTEGER ,       INTENT(OUT) :: Status
3065   INTEGER          :: idx
3066   INTEGER          :: stat
3067   CHARACTER*(1000) :: VALUE
3069   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_td_integer')
3071   Status = WRF_NO_ERR
3072   
3073   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), DateStr, &
3074        "none", Value,stat)
3075   if (stat /= 0) then
3076      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
3077      Status = WRF_WARN_VAR_NF
3078      RETURN
3079   endif
3081   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
3082   if (stat .ne. 0) then
3083      CALL wrf_message("Reading data from"//Value//"failed")
3084      Status = WRF_WARN_COUNT_TOO_LONG
3085      RETURN
3086   endif
3087   Outcount = idx
3089   RETURN
3090 END SUBROUTINE ext_gr1_get_dom_td_integer 
3092 !*****************************************************************************
3094 SUBROUTINE ext_gr1_get_dom_td_logical ( DataHandle,Element, DateStr,  Data, &
3095      Count, Outcount, Status )
3097   USE gr1_data_info
3098   IMPLICIT NONE
3099 #include "wrf_status_codes.h"
3100   INTEGER ,       INTENT(IN)  :: DataHandle
3101   CHARACTER*(*) :: Element
3102   CHARACTER*(*) :: DateStr
3103   logical ,       INTENT(OUT) :: Data(*)
3104   INTEGER ,       INTENT(IN)  :: Count
3105   INTEGER ,       INTENT(OUT) :: OutCount
3106   INTEGER ,       INTENT(OUT) :: Status
3107   INTEGER          :: idx
3108   INTEGER          :: stat
3109   CHARACTER*(1000) :: VALUE
3111   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_td_logical')
3113   Status = WRF_NO_ERR
3114   
3115   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), DateStr, &
3116        "none", Value, stat)
3117   if (stat /= 0) then
3118      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
3119      Status = WRF_WARN_VAR_NF
3120      RETURN
3121   endif
3123   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
3124   if (stat .ne. 0) then
3125      CALL wrf_message("Reading data from"//Value//"failed")
3126      Status = WRF_WARN_COUNT_TOO_LONG
3127      RETURN
3128   endif
3129   Outcount = idx
3131   RETURN
3132 END SUBROUTINE ext_gr1_get_dom_td_logical 
3134 !*****************************************************************************
3136 SUBROUTINE ext_gr1_get_dom_td_char ( DataHandle,Element, DateStr,  Data,  &
3137      Status )
3139   USE gr1_data_info
3140   IMPLICIT NONE
3141 #include "wrf_status_codes.h"
3142   INTEGER ,       INTENT(IN)  :: DataHandle
3143   CHARACTER*(*) :: Element
3144   CHARACTER*(*) :: DateStr
3145   CHARACTER*(*) :: Data
3146   INTEGER ,       INTENT(OUT) :: Status
3147   INTEGER       :: stat
3149   Status = WRF_NO_ERR
3150   
3151   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_td_char')
3153   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), DateStr, &
3154        "none", Data, stat)
3155   if (stat /= 0) then
3156      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
3157      Status = WRF_WARN_VAR_NF
3158      RETURN
3159   endif
3161   RETURN
3162 END SUBROUTINE ext_gr1_get_dom_td_char 
3164 !*****************************************************************************
3166 SUBROUTINE ext_gr1_get_dom_td_double ( DataHandle,Element, DateStr,  Data, &
3167      Count, Outcount, Status )
3168   USE gr1_data_info
3169   IMPLICIT NONE
3170 #include "wrf_status_codes.h"
3171   INTEGER ,       INTENT(IN)  :: DataHandle
3172   CHARACTER*(*) , INTENT(IN)  :: Element
3173   CHARACTER*(*) , INTENT(IN)  :: DateStr
3174   real*8 ,            INTENT(OUT) :: Data(*)
3175   INTEGER ,       INTENT(IN)  :: Count
3176   INTEGER ,       INTENT(OUT)  :: OutCount
3177   INTEGER ,       INTENT(OUT) :: Status
3178   INTEGER          :: idx
3179   INTEGER          :: stat
3180   CHARACTER*(1000) :: VALUE
3182   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_td_double')
3184   Status = WRF_NO_ERR
3185   
3186   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), DateStr, &
3187        "none", Value, stat)
3188   if (stat /= 0) then
3189      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
3190      Status = WRF_WARN_VAR_NF
3191      RETURN
3192   endif
3194   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
3195   if (stat .ne. 0) then
3196      CALL wrf_message("Reading data from"//Value//"failed")
3197      Status = WRF_WARN_COUNT_TOO_LONG
3198      RETURN
3199   endif
3200   Outcount = idx
3202 RETURN
3203 END SUBROUTINE ext_gr1_get_dom_td_double
3205 !******************************************************************************
3206 !* End of get_dom_td_* routines
3207 !******************************************************************************
3210 !******************************************************************************
3211 !* Start of put_dom_td_* routines
3212 !******************************************************************************
3215 SUBROUTINE ext_gr1_put_dom_td_real8 ( DataHandle,Element, DateStr,  Data, &
3216      Count,  Status )
3218   USE gr1_data_info
3219   IMPLICIT NONE
3220 #include "wrf_status_codes.h"
3221   INTEGER ,       INTENT(IN)  :: DataHandle
3222   CHARACTER*(*) :: Element
3223   CHARACTER*(*) :: DateStr
3224   real*8 ,        INTENT(IN)  :: Data(*)
3225   INTEGER ,       INTENT(IN)  :: Count
3226   INTEGER ,       INTENT(OUT) :: Status
3227   CHARACTER(len=1000) :: tmpstr(1000)
3228   INTEGER             :: idx
3230   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_td_real8')
3232   if (committed(DataHandle)) then
3234      do idx = 1,Count
3235         write(tmpstr(idx),'(G17.10)')Data(idx)
3236      enddo
3238      CALL gr1_build_string (td_output(DataHandle), DateStr//';'//Element, tmpstr, &
3239           Count, Status)
3241   endif
3243   RETURN
3244 END SUBROUTINE ext_gr1_put_dom_td_real8 
3246 !*****************************************************************************
3248 SUBROUTINE ext_gr1_put_dom_td_integer ( DataHandle,Element, DateStr,  Data, &
3249      Count,  Status )
3251   USE gr1_data_info
3252   IMPLICIT NONE
3253 #include "wrf_status_codes.h"
3254   INTEGER ,       INTENT(IN)  :: DataHandle
3255   CHARACTER*(*) :: Element
3256   CHARACTER*(*) :: DateStr
3257   integer ,       INTENT(IN)  :: Data(*)
3258   INTEGER ,       INTENT(IN)  :: Count
3259   INTEGER ,       INTENT(OUT) :: Status
3260   CHARACTER(len=1000) :: tmpstr(1000)
3261   INTEGER             :: idx
3263   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_td_integer')
3265   if (committed(DataHandle)) then
3267      do idx = 1,Count
3268         write(tmpstr(idx),'(G17.10)')Data(idx)
3269      enddo
3270      
3271      CALL gr1_build_string (td_output(DataHandle), DateStr//';'//Element, tmpstr, &
3272           Count, Status)
3274   endif
3276   RETURN
3277 END SUBROUTINE ext_gr1_put_dom_td_integer
3279 !*****************************************************************************
3281 SUBROUTINE ext_gr1_put_dom_td_logical ( DataHandle,Element, DateStr,  Data, &
3282      Count,  Status )
3284   USE gr1_data_info
3285   IMPLICIT NONE
3286 #include "wrf_status_codes.h"
3287   INTEGER ,       INTENT(IN)  :: DataHandle
3288   CHARACTER*(*) :: Element
3289   CHARACTER*(*) :: DateStr
3290   logical ,       INTENT(IN)  :: Data(*)
3291   INTEGER ,       INTENT(IN)  :: Count
3292   INTEGER ,       INTENT(OUT) :: Status
3293   CHARACTER(len=1000) :: tmpstr(1000)
3294   INTEGER             :: idx
3296   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_td_logical')
3298   if (committed(DataHandle)) then
3300      do idx = 1,Count
3301         write(tmpstr(idx),'(G17.10)')Data(idx)
3302      enddo
3303      
3304      CALL gr1_build_string (td_output(DataHandle), DateStr//';'//Element, tmpstr, &
3305           Count, Status)
3307   endif
3309   RETURN
3310 END SUBROUTINE ext_gr1_put_dom_td_logical
3312 !*****************************************************************************
3314 SUBROUTINE ext_gr1_put_dom_td_char ( DataHandle,Element, DateStr,  Data, &
3315      Status )
3317   USE gr1_data_info
3318   IMPLICIT NONE
3319 #include "wrf_status_codes.h"
3320   INTEGER ,       INTENT(IN)  :: DataHandle
3321   CHARACTER*(*) :: Element
3322   CHARACTER*(*) :: DateStr
3323   CHARACTER(len=*), INTENT(IN)  :: Data
3324   INTEGER ,       INTENT(OUT) :: Status
3325   CHARACTER(len=1000) :: tmpstr(1)
3327   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_td_char')
3329   if (committed(DataHandle)) then
3331      write(tmpstr(1),*)Data
3333      CALL gr1_build_string (td_output(DataHandle), DateStr//';'//Element, tmpstr, &
3334           1, Status)
3336   endif
3338   RETURN
3339 END SUBROUTINE ext_gr1_put_dom_td_char 
3341 !*****************************************************************************
3343 SUBROUTINE ext_gr1_put_dom_td_double ( DataHandle,Element, DateStr,  Data, &
3344      Count,  Status )
3345   USE gr1_data_info
3346   IMPLICIT NONE
3347 #include "wrf_status_codes.h"
3348   INTEGER ,       INTENT(IN)  :: DataHandle
3349   CHARACTER*(*) , INTENT(IN)  :: Element
3350   CHARACTER*(*) , INTENT(IN)  :: DateStr
3351   real*8 ,            INTENT(IN) :: Data(*)
3352   INTEGER ,       INTENT(IN)  :: Count
3353   INTEGER ,       INTENT(OUT) :: Status
3354   CHARACTER(len=1000) :: tmpstr(1000)
3355   INTEGER             :: idx
3357   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_td_double')
3359   if (committed(DataHandle)) then
3361      do idx = 1,Count
3362         write(tmpstr(idx),'(G17.10)')Data(idx)
3363      enddo
3365      CALL gr1_build_string (td_output(DataHandle), DateStr//';'//Element, tmpstr, &
3366           Count, Status)
3368   endif
3370 RETURN
3371 END SUBROUTINE ext_gr1_put_dom_td_double
3373 !*****************************************************************************
3375 SUBROUTINE ext_gr1_put_dom_td_real ( DataHandle,Element, DateStr,  Data, &
3376      Count,  Status )
3378   USE gr1_data_info
3379   IMPLICIT NONE
3380 #include "wrf_status_codes.h"
3381   INTEGER ,       INTENT(IN)  :: DataHandle
3382   CHARACTER*(*) :: Element
3383   CHARACTER*(*) :: DateStr
3384   real ,          INTENT(IN)  :: Data(*)
3385   INTEGER ,       INTENT(IN)  :: Count
3386   INTEGER ,       INTENT(OUT) :: Status
3387   CHARACTER(len=1000) :: tmpstr(1000)
3388   INTEGER             :: idx
3390   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_td_real')
3392   if (committed(DataHandle)) then
3394      do idx = 1,Count
3395         write(tmpstr(idx),'(G17.10)')Data(idx)
3396      enddo
3397      
3398      CALL gr1_build_string (td_output(DataHandle), DateStr//';'//Element, tmpstr, &
3399           Count, Status)
3401   endif
3403   RETURN
3404 END SUBROUTINE ext_gr1_put_dom_td_real 
3407 !******************************************************************************
3408 !* End of put_dom_td_* routines
3409 !******************************************************************************
3412 !*****************************************************************************
3414 SUBROUTINE gr1_build_string (string, Element, Value, Count, Status)
3416   IMPLICIT NONE
3417 #include "wrf_status_codes.h"
3419   CHARACTER (LEN=*) , INTENT(INOUT) :: string
3420   CHARACTER (LEN=*) , INTENT(IN)    :: Element
3421   CHARACTER (LEN=*) , INTENT(IN)    :: Value(*)
3422   INTEGER ,           INTENT(IN)    :: Count
3423   INTEGER ,           INTENT(OUT)   :: Status
3425   CHARACTER (LEN=2)                 :: lf
3426   INTEGER                           :: IDX
3428   lf=char(10)//' '
3429   if (len_trim(string) == 0) then
3430      string = lf//Element//' = '
3431   else
3432      string = trim(string)//lf//Element//' = '
3433   endif
3434   do idx = 1,Count
3435      if (idx > 1) then
3436         string = trim(string)//','
3437      endif
3438      string = trim(string)//' '//trim(adjustl(Value(idx)))
3439   enddo
3441   Status = WRF_NO_ERR
3443 END SUBROUTINE gr1_build_string
3445 !*****************************************************************************
3447 SUBROUTINE gr1_get_new_handle(DataHandle)
3448   USE gr1_data_info
3449   IMPLICIT NONE
3450   
3451   INTEGER ,       INTENT(OUT)  :: DataHandle
3452   INTEGER :: i
3454   DataHandle = -1
3455   do i=firstFileHandle, maxFileHandles
3456      if (.NOT. used(i)) then
3457         DataHandle = i
3458         used(i) = .true.
3459         exit
3460      endif
3461   enddo
3463   RETURN
3464 END SUBROUTINE gr1_get_new_handle
3467 !******************************************************************************
3470 SUBROUTINE gr1_get_levels(VarName, zidx, zsize, soil_layers, vert_stag, fraction, &
3471      is_press_levels, is_turb_layers, vert_unit, level1, level2)
3473   use gr1_data_info
3474   IMPLICIT NONE
3476   integer :: zidx
3477   integer :: zsize
3478   logical :: soil_layers
3479   logical :: vert_stag
3480   logical :: fraction
3481   logical :: is_press_levels
3482   logical :: is_turb_layers 
3483   integer :: vert_unit
3484   integer :: level1
3485   integer :: level2
3486   character (LEN=*) :: VarName
3488   ! Setup vert_unit, and vertical levels in grib units
3490   if ((VarName .eq. 'LANDUSEF') .or. (VarName .eq. 'SOILCTOP') &
3491        .or. (VarName .eq. 'SOILCBOT')) then
3492      vert_unit = 109;
3493      level1 = zidx
3494      level2 = 0
3495   else if ((zsize .gt. 1) .and. (.not. soil_layers) .and. (.not. fraction) &
3496             .and. (.not. is_press_levels) .and. (.not. is_turb_layers)) then
3497      vert_unit = 119;
3498      if (vert_stag) then
3499         level1 = (10000*full_eta(zidx)+0.5)
3500      else
3501         level1 = (10000*half_eta(zidx)+0.5)
3502      endif
3503      level2 = 0
3504   else
3505      ! Set the vertical coordinate and level for soil and 2D fields
3506      if (fraction) then
3507         vert_unit = 109
3508         level1 = zidx
3509         level2 = 0           
3510      else if (soil_layers) then
3511         vert_unit = 112
3512         level1 = 100*(soil_depth(zidx) - 0.5*soil_thickness(zidx))+0.5
3513         level2 = 100*(soil_depth(zidx) + 0.5*soil_thickness(zidx))+0.5
3514      ! Added pressure level as vertical unit - GAC 20140402
3515      else if (is_press_levels) then
3516         vert_unit = 100
3517         level1 = press_levels(zidx)/100.
3518         level2 = 0
3519      else if (is_turb_layers) then
3520         vert_unit = 106
3521         level1 = turb_layer_top(zidx)/100.
3522         level2 = turb_layer_bot(zidx)/100.
3523      else if (VarName .eq. 'mu') then
3524         vert_unit = 200
3525         level1 = 0
3526         level2 = 0
3527      else if ((VarName .eq. 'Q2') .or. (VarName .eq. 'TH2') .or. &
3528         (VarName .eq. 'T2')) then
3529         vert_unit = 105
3530         level1 = 2
3531         level2 = 0
3532      else if ((VarName .eq. 'Q10') .or. (VarName .eq. 'TH10') .or. &
3533           (VarName .eq. 'U10') .or. (VarName .eq. 'V10')) then
3534         vert_unit = 105
3535         level1 = 10
3536         level2 = 0
3537      ! Until a more sophisticated way to do this becomes clear,
3538      ! adding AFWA diagnostic variable vertical unit and level info
3539      ! here.  GAC 20140402
3540      else if ((VarName .eq. 'TCOLI_MAX') .or.                      &
3541              (Varname .eq. 'GRPL_FLX_MAX') .or.                    &
3542              (Varname .eq. 'VIL') .or.                             &
3543              (Varname .eq. 'RADARVIL') .or.                        &
3544              (Varname .eq. 'FZLEV') .or.                           &
3545              (Varname .eq. 'REFD_COM') .or.                        &
3546              (Varname .eq. 'ICINGTOP') .or.                        &
3547              (Varname .eq. 'ICINGBOT') .or.                        &
3548              (Varname .eq. 'ICING_LG') .or.                        &
3549              (Varname .eq. 'ICING_SM') .or.                        &
3550              (Varname .eq. 'AFWA_CLOUD') .or.                      &
3551              (Varname .eq. 'AFWA_CLOUD_CEIL') .or.                 &
3552              (Varname .eq. 'AFWA_HAIL') .or.                       &
3553              (Varname .eq. 'AFWA_TORNADO') .or.                    &
3554              (Varname .eq. 'AFWA_PWAT') .or.                       &
3555              (Varname .eq. 'QICING_LG_MAX') .or.                   &
3556              (Varname .eq. 'QICING_SM_MAX')) then
3557         vert_unit = 200
3558         level1 = 0
3559         level2 = 0
3560      else if (VarName .eq. 'REFD_MAX') then
3561         vert_unit = 105
3562         level1 = 1000
3563         level2 = 0
3564      else if (VarName .eq. 'WSPD10MAX') then
3565         vert_unit = 105
3566         level1 = 10
3567         level2 = 0
3568      else if (VarName .eq. 'UP_HELI_MAX') then
3569         vert_unit = 106
3570         level1 = 50 !5000 m
3571         level2 = 20 !2000 m
3572      else if (Varname .eq. 'AFWA_LLWS') then
3573         vert_unit = 106
3574         level1 = 20 !2000 m
3575         level2 = 0  !0 m
3576      else if ((Varname .eq. 'AFWA_LLTURB') .or.                   &
3577               (Varname .eq. 'AFWA_LLTURBLGT') .or.                &
3578               (Varname .eq. 'AFWA_LLTURBMDT') .or.                &
3579               (Varname .eq. 'AFWA_LLTURBSVR')) then
3580         vert_unit = 106
3581         level1 = 15 !1500 m
3582         level2 = 0  !0 m
3583      else if ((VarName .eq. 'W_UP_MAX') .or. (VarName .eq. 'W_DN_MAX')) then
3584         vert_unit = 101
3585         level1 = 40  !400 mb
3586         level2 = 100 !1000 mb
3587      else if (VarName .eq. 'AFWA_LIDX') then
3588         vert_unit = 101
3589         level1 = 50  !500 mb
3590         level2 = 100 !1000 mb
3591      else if (VarName .eq. 'AFWA_MSLP') then
3592         vert_unit = 102  ! Mean sea level
3593         level1 = 0
3594         level2 = 0
3595      else if ((VarName .eq. 'AFWA_CAPE_MU') .or.                  &
3596               (VarName .eq. 'AFWA_CIN_MU')) then
3597         vert_unit = 116  ! Pressure above ground
3598         !level1 = 180     ! 180 mb AGL
3599         !level2 = 0
3600         level1 = 0       ! 180 mb AGL
3601         level2 = 180
3602      else 
3603         vert_unit = 1
3604         level1 = 0
3605         level2 = 0
3606      endif
3607   endif
3609 end SUBROUTINE gr1_get_levels
3611 !*****************************************************************************
3614 SUBROUTINE gr1_fill_eta_levels(fileindex, FileFd, grib_tables, VarName, eta_levels)
3615   IMPLICIT NONE
3617   CHARACTER (len=*) :: fileindex
3618   INTEGER   :: FileFd
3619   CHARACTER (len=*) :: grib_tables
3620   character (len=*) :: VarName
3621   REAL,DIMENSION(*) :: eta_levels
3623   INTEGER   :: center, subcenter, parmtbl
3624   INTEGER   :: swapped
3625   INTEGER   :: leveltype
3626   INTEGER   :: idx
3627   INTEGER   :: parmid
3628   INTEGER   :: tablenum
3629   REAL      :: tmp
3630   INTEGER   :: numindices
3631   integer , DIMENSION(1000)   :: indices
3633   !
3634   ! Read the levels from the grib file
3635   !
3636   CALL GET_GRIB_PARAM(grib_tables, VarName, center, subcenter, parmtbl, &
3637        tablenum, parmid)
3639   if (parmid == -1) then
3640      call wrf_message ('Error getting grib parameter')
3641   endif
3643   leveltype = 119
3645   CALL GET_GRIB_INDICES(fileindex(:), center, subcenter, parmtbl, &
3646        parmid, "*", leveltype, &
3647        -HUGE(1), -HUGE(1), -HUGE(1), -HUGE(1), indices, numindices)
3649   do idx = 1,numindices
3650      CALL READ_GRIB(fileindex(:),FileFd,indices(idx),eta_levels(idx))
3651   enddo
3653   !
3654   ! Sort the levels--from highest (bottom) to lowest (top)
3655   !
3656   swapped = 1
3657   sortloop : do
3658      if (swapped /= 1) exit sortloop
3659      swapped = 0
3660      do idx=2, numindices
3661         !
3662         ! Remove duplicate levels, caused by multiple time periods in a 
3663         ! single file.
3664         !
3665         if (eta_levels(idx) == eta_levels(idx-1)) eta_levels(idx) = 0.0
3666         if (eta_levels(idx) > eta_levels(idx-1)) then
3667           tmp = eta_levels(idx)
3668           eta_levels(idx) = eta_levels(idx - 1)
3669           eta_levels(idx - 1) = tmp
3670           swapped = 1
3671         endif
3672      enddo
3673   enddo sortloop
3675 end subroutine gr1_fill_eta_levels