standard WRF version 3.0.1.1
[wrffire.git] / wrfv2_fire / external / io_grib1 / io_grib1.F
blobcdf5badd9501b6248fbb62ce0eb964ddf0d99b4f
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   = 200
31   integer                , parameter       :: maxLevels        = 1000
32   integer                , parameter       :: maxSoilLevels    = 100
33   integer                , parameter       :: maxDomains       = 500
35   logical ,      dimension(maxFileHandles) :: committed, opened, used
36   character*128, dimension(maxFileHandles) :: DataFile
37   integer,       dimension(maxFileHandles) :: FileFd
38   integer,       dimension(maxFileHandles) :: FileStatus
39   REAL,          dimension(maxLevels)      :: half_eta, full_eta
40   REAL,          dimension(maxSoilLevels)  :: soil_depth, soil_thickness
41   character*24                             :: StartDate = ''
42   character*24                             :: InputProgramName = ''
43   integer                                  :: projection
44   integer                                  :: wg_grid_id
45   real                                     :: dx,dy
46   real                                     :: truelat1, truelat2
47   real                                     :: center_lat, center_lon
48   real                                     :: proj_central_lon
49   real                                     :: timestep
50   character,     dimension(:), pointer     :: grib_tables
51   logical                                  :: table_filled = .FALSE.
52   character,     dimension(:), pointer     :: grid_info
53   integer                                  :: full_xsize, full_ysize
54   integer, dimension(maxDomains)           :: domains = -1
55   integer                                  :: max_domain = 0
56   
57   TYPE :: HandleVar
58      character, dimension(:), pointer      :: fileindex(:)
59      integer                               :: CurrentTime
60      integer                               :: NumberTimes
61      character (DateStrLen), dimension(:),pointer  :: Times(:)
62   ENDTYPE
63   TYPE (HandleVar), dimension(maxFileHandles) :: fileinfo
65   TYPE :: prevdata
66      integer :: fcst_secs_rainc
67      integer :: fcst_secs_rainnc
68      real, dimension(:,:), pointer         :: rainc, rainnc
69   END TYPE prevdata
70   TYPE (prevdata), DIMENSION(500)          :: lastdata
72   TYPE :: initdata
73      real,         dimension(:,:), pointer :: snod
74   END TYPE initdata
76   TYPE (initdata), dimension(maxDomains)   :: firstdata
78   TYPE :: prestype
79      real,         dimension(:,:,:), pointer :: vals
80      logical                                :: newtime
81      character*120                          :: lastDateStr
82   END TYPE prestype
84   character*120, dimension(maxDomains)     :: lastDateStr
86   TYPE (prestype), dimension(maxDomains)   :: pressure
87   TYPE (prestype), dimension(maxDomains)   :: geopotential
89   integer                                  :: center, subcenter, parmtbl
91   character(len=30000), dimension(maxFileHandles) :: td_output
92   character(len=30000), dimension(maxFileHandles) :: ti_output
94   logical                                  :: WrfIOnotInitialized = .true.
96 end module gr1_data_info
99 subroutine ext_gr1_ioinit(SysDepInfo,Status)
101   USE gr1_data_info
102   implicit none
103 #include "wrf_status_codes.h"
104 #include "wrf_io_flags.h"
105   CHARACTER*(*), INTENT(IN) :: SysDepInfo
106   integer ,intent(out) :: Status
107   integer :: i
108   integer :: size, istat
109   CHARACTER (LEN=300) :: wrf_err_message
111   call wrf_debug ( DEBUG , 'Entering ext_gr1_ioinit')
113   do i=firstFileHandle, maxFileHandles
114         used(i) = .false.
115         committed(i) = .false.
116         opened(i) = .false.
117         td_output(i) = ''
118         ti_output(i) = ''
119   enddo
120   domains(:) = -1
122   do i = 1, maxDomains
123     pressure(i)%newtime = .false.
124     pressure(i)%lastDateStr = ''
125     geopotential(i)%newtime = .false.
126     geopotential(i)%lastDateStr = ''
127     lastDateStr(i) = ''
128   enddo
130   lastdata%fcst_secs_rainc = 0
131   lastdata%fcst_secs_rainnc = 0
132   FileStatus(1:maxFileHandles) = WRF_FILE_NOT_OPENED
133   WrfIOnotInitialized = .false.
135   Status = WRF_NO_ERR
137   return
138 end subroutine ext_gr1_ioinit
140 !*****************************************************************************
142 subroutine ext_gr1_ioexit(Status)
144   USE gr1_data_info
145   implicit none
146 #include "wrf_status_codes.h"
147   integer istat
148   integer ,intent(out) :: Status
150   call wrf_debug ( DEBUG , 'Entering ext_gr1_ioexit')
152   if (table_filled) then
153      CALL free_gribmap(grib_tables)
154      DEALLOCATE(grib_tables, stat=istat)
155      table_filled = .FALSE.
156   endif
157   IF ( ASSOCIATED ( grid_info ) ) THEN
158     DEALLOCATE(grid_info, stat=istat)
159   ENDIF
160   NULLIFY(grid_info)
162   Status = WRF_NO_ERR
164   return
165 end subroutine ext_gr1_ioexit
167 !*****************************************************************************
169 SUBROUTINE ext_gr1_open_for_read_begin ( FileName , Comm_compute, Comm_io, &
170      SysDepInfo, DataHandle , Status )
172   USE gr1_data_info
173   IMPLICIT NONE
174 #include "wrf_status_codes.h"
175 #include "wrf_io_flags.h"
176   CHARACTER*(*) :: FileName
177   INTEGER ,       INTENT(IN)  :: Comm_compute , Comm_io
178   CHARACTER*(*) :: SysDepInfo
179   INTEGER ,       INTENT(OUT) :: DataHandle
180   INTEGER ,       INTENT(OUT) :: Status
181   integer                     :: ierr
182   integer                     :: size
183   integer                     :: idx
184   integer                     :: parmid
185   integer                     :: dpth_parmid
186   integer                     :: thk_parmid
187   integer                     :: leveltype
188   integer , DIMENSION(1000)   :: indices
189   integer                     :: numindices
190   real , DIMENSION(1000)      :: levels
191   real                        :: tmp
192   integer                     :: swapped
193   integer                     :: etaidx
194   integer                     :: grb_index
195   integer                     :: level1, level2
196   integer   :: tablenum
197   integer   :: stat
198   integer   :: endchar
199   integer   :: last_grb_index
200   CHARACTER (LEN=300) :: wrf_err_message
202   call wrf_debug ( DEBUG , 'Entering ext_gr1_open_for_read_begin')
204   CALL gr1_get_new_handle(DataHandle)
206   if (DataHandle .GT. 0) then
207      CALL open_file(TRIM(FileName), 'r', FileFd(DataHandle), ierr)
208      if (ierr .ne. 0) then
209         Status = WRF_ERR_FATAL_BAD_FILE_STATUS
210      else
211         opened(DataHandle) = .true.
212         DataFile(DataHandle) = TRIM(FileName)
213         FileStatus(DataHandle) = WRF_FILE_OPENED_NOT_COMMITTED
214      endif
215   else
216      Status = WRF_WARN_TOO_MANY_FILES
217      return
218   endif
220   ! Read the grib index file first
221   if (.NOT. table_filled) then
222      table_filled = .TRUE.
223      CALL GET_GRIB1_TABLES_SIZE(size)
224      ALLOCATE(grib_tables(1:size), STAT=ierr)
225      CALL LOAD_GRIB1_TABLES ("gribmap.txt", grib_tables, ierr)
226      if (ierr .ne. 0) then
227         DEALLOCATE(grib_tables)
228         WRITE( wrf_err_message , * ) &
229              'Could not open file gribmap.txt '
230         CALL wrf_error_fatal ( TRIM ( wrf_err_message ) )
231         Status = WRF_ERR_FATAL_BAD_FILE_STATUS
232         return
233      endif
234   endif
236   ! Begin by indexing file and reading metadata into structure.
237   CALL GET_FILEINDEX_SIZE(size)
238   ALLOCATE(fileinfo(DataHandle)%fileindex(1:size), STAT=ierr)
240   CALL ALLOC_INDEX_FILE(fileinfo(DataHandle)%fileindex(:))
241   CALL INDEX_FILE(FileFd(DataHandle),fileinfo(DataHandle)%fileindex(:))
243   ! Get times into Times variable
244   CALL GET_NUM_TIMES(fileinfo(DataHandle)%fileindex(:), &
245        fileinfo(DataHandle)%NumberTimes);
247   ALLOCATE(fileinfo(DataHandle)%Times(1:fileinfo(DataHandle)%NumberTimes), STAT=ierr)
248   do idx = 1,fileinfo(DataHandle)%NumberTimes
249      CALL GET_TIME(fileinfo(DataHandle)%fileindex(:),idx, &
250           fileinfo(DataHandle)%Times(idx))
251   enddo
253   ! CurrentTime starts as 0.  The first time in the file is 1.  So,
254   !   until set_time or get_next_time is called, the current time
255   !   is not set.
256   fileinfo(DataHandle)%CurrentTime = 0
258   CALL gr1_fill_eta_levels(fileinfo(DataHandle)%fileindex(:), &
259        FileFd(DataHandle), & 
260        grib_tables, "ZNW", full_eta)
261   CALL gr1_fill_eta_levels(fileinfo(DataHandle)%fileindex(:), FileFd(DataHandle), &
262        grib_tables, "ZNU", half_eta)
264   ! 
265   ! Now, get the soil levels
266   !
267   CALL GET_GRIB_PARAM(grib_tables, "ZS", center, subcenter, parmtbl, &
268        tablenum, dpth_parmid)
269   CALL GET_GRIB_PARAM(grib_tables,"DZS", center, subcenter, parmtbl, &
270        tablenum, thk_parmid)
271   if (dpth_parmid == -1) then
272      call wrf_message ('Error getting grib parameter')
273   endif
275   leveltype = 112
277   CALL GET_GRIB_INDICES(fileinfo(DataHandle)%fileindex(:),center, subcenter, parmtbl, &
278        dpth_parmid,"*",leveltype, &
279        -HUGE(1),-HUGE(1), -HUGE(1),-HUGE(1),indices,numindices)
281   last_grb_index = -1;
282   do idx = 1,numindices
283      CALL READ_GRIB(fileinfo(DataHandle)%fileindex(:), FileFd(DataHandle), &
284           indices(idx), soil_depth(idx))
285      !
286      ! Now read the soil thickenesses
287      !
288      CALL GET_LEVEL1(fileinfo(DataHandle)%fileindex(:),indices(idx),level1)
289      CALL GET_LEVEL2(fileinfo(DataHandle)%fileindex(:),indices(idx),level2)
290      CALL GET_GRIB_INDEX_GUESS(fileinfo(DataHandle)%fileindex(:), &
291           center, subcenter, parmtbl, thk_parmid,"*",leveltype, &
292           level1,level2,-HUGE(1),-HUGE(1), last_grb_index+1, grb_index)
293      CALL READ_GRIB(fileinfo(DataHandle)%fileindex(:),FileFd(DataHandle),grb_index, &
294           soil_thickness(idx))
296      last_grb_index = grb_index
297   enddo
298   
301   !
302   ! Fill up any variables that need to be retrieved from Metadata
303   !
304   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), 'PROGRAM_NAME', "none", &
305        "none", InputProgramName, stat)
306   if (stat /= 0) then
307      CALL wrf_debug (DEBUG , "PROGRAM_NAME not found in input METADATA")
308   else 
309      endchar = SCAN(InputProgramName," ")
310      InputProgramName = InputProgramName(1:endchar)
311   endif
313   call wrf_debug ( DEBUG , 'Exiting ext_gr1_open_for_read_begin')
315   RETURN
316 END SUBROUTINE ext_gr1_open_for_read_begin
318 !*****************************************************************************
320 SUBROUTINE ext_gr1_open_for_read_commit( DataHandle , Status )
322   USE gr1_data_info
323   IMPLICIT NONE
324 #include "wrf_status_codes.h"
325 #include "wrf_io_flags.h"
326   character(len=1000) :: msg
327   INTEGER ,       INTENT(IN ) :: DataHandle
328   INTEGER ,       INTENT(OUT) :: Status
330   call wrf_debug ( DEBUG , 'Entering ext_gr1_open_for_read_commit')
332   Status = WRF_NO_ERR
333   if(WrfIOnotInitialized) then
334     Status = WRF_IO_NOT_INITIALIZED
335     write(msg,*) 'ext_gr1_ioinit was not called ',__FILE__,', line', __LINE__
336     call wrf_debug ( FATAL , msg)
337     return
338   endif
339   committed(DataHandle) = .true.
340   FileStatus(DataHandle) = WRF_FILE_OPENED_FOR_READ
342   Status = WRF_NO_ERR
344   RETURN
345 END SUBROUTINE ext_gr1_open_for_read_commit
347 !*****************************************************************************
349 SUBROUTINE ext_gr1_open_for_read ( FileName , Comm_compute, Comm_io, &
350      SysDepInfo, DataHandle , Status )
352   USE gr1_data_info
353   IMPLICIT NONE
354 #include "wrf_status_codes.h"
355 #include "wrf_io_flags.h"
356   CHARACTER*(*) :: FileName
357   INTEGER ,       INTENT(IN)  :: Comm_compute , Comm_io
358   CHARACTER*(*) :: SysDepInfo
359   INTEGER ,       INTENT(OUT) :: DataHandle
360   INTEGER ,       INTENT(OUT) :: Status
363   call wrf_debug ( DEBUG , 'Entering ext_gr1_open_for_read')
365   DataHandle = 0   ! dummy setting to quiet warning message
366   CALL ext_gr1_open_for_read_begin( FileName, Comm_compute, Comm_io, &
367        SysDepInfo, DataHandle, Status )
368   IF ( Status .EQ. WRF_NO_ERR ) THEN
369      FileStatus(DataHandle) = WRF_FILE_OPENED_NOT_COMMITTED
370      CALL ext_gr1_open_for_read_commit( DataHandle, Status )
371   ENDIF
372   return
374   RETURN  
375 END SUBROUTINE ext_gr1_open_for_read
377 !*****************************************************************************
379 SUBROUTINE ext_gr1_open_for_write_begin(FileName, Comm, IOComm, SysDepInfo, &
380      DataHandle, Status)
381   
382   USE gr1_data_info
383   implicit none
384 #include "wrf_status_codes.h"
385 #include "wrf_io_flags.h"
387   character*(*)        ,intent(in)  :: FileName
388   integer              ,intent(in)  :: Comm
389   integer              ,intent(in)  :: IOComm
390   character*(*)        ,intent(in)  :: SysDepInfo
391   integer              ,intent(out) :: DataHandle
392   integer              ,intent(out) :: Status
393   integer :: ierr
394   CHARACTER (LEN=300) :: wrf_err_message
395   integer             :: size
397   call wrf_debug ( DEBUG , 'Entering ext_gr1_open_for_write_begin')
399   if (.NOT. table_filled) then
400      table_filled = .TRUE.
401      CALL GET_GRIB1_TABLES_SIZE(size)
402      ALLOCATE(grib_tables(1:size), STAT=ierr)
403      CALL LOAD_GRIB1_TABLES ("gribmap.txt", grib_tables, ierr)
404      if (ierr .ne. 0) then
405         DEALLOCATE(grib_tables)
406         WRITE( wrf_err_message , * ) &
407              'Could not open file gribmap.txt '
408         CALL wrf_error_fatal ( TRIM ( wrf_err_message ) )
409         Status = WRF_ERR_FATAL_BAD_FILE_STATUS
410         return
411      endif
412   endif
414   Status = WRF_NO_ERR
415   CALL gr1_get_new_handle(DataHandle)
416   if (DataHandle .GT. 0) then
417      CALL open_file(TRIM(FileName), 'w', FileFd(DataHandle), ierr)
418      if (ierr .ne. 0) then
419         Status = WRF_WARN_WRITE_RONLY_FILE
420      else
421         opened(DataHandle) = .true.
422         DataFile(DataHandle) = TRIM(FileName)
423         FileStatus(DataHandle) = WRF_FILE_OPENED_NOT_COMMITTED
424      endif
425      committed(DataHandle) = .false.
426      td_output(DataHandle) = ''
427   else
428      Status = WRF_WARN_TOO_MANY_FILES
429   endif
431   RETURN  
432 END SUBROUTINE ext_gr1_open_for_write_begin
434 !*****************************************************************************
436 SUBROUTINE ext_gr1_open_for_write_commit( DataHandle , Status )
438   USE gr1_data_info
439   IMPLICIT NONE
440 #include "wrf_status_codes.h"
441 #include "wrf_io_flags.h"
442   INTEGER ,       INTENT(IN ) :: DataHandle
443   INTEGER ,       INTENT(OUT) :: Status
445   call wrf_debug ( DEBUG , 'Entering ext_gr1_open_for_write_commit')
447   IF ( opened( DataHandle ) ) THEN
448     IF ( used( DataHandle ) ) THEN
449       committed(DataHandle) = .true.
450       FileStatus(DataHandle) = WRF_FILE_OPENED_FOR_WRITE
451     ENDIF
452   ENDIF
454   Status = WRF_NO_ERR
456   RETURN  
457 END SUBROUTINE ext_gr1_open_for_write_commit
459 !*****************************************************************************
461 subroutine ext_gr1_inquiry (Inquiry, Result, Status)
462   use gr1_data_info
463   implicit none
464 #include "wrf_status_codes.h"
465   character *(*), INTENT(IN)    :: Inquiry
466   character *(*), INTENT(OUT)   :: Result
467   integer        ,INTENT(INOUT) :: Status
468   SELECT CASE (Inquiry)
469   CASE ("RANDOM_WRITE","RANDOM_READ")
470      Result='ALLOW'
471   CASE ("SEQUENTIAL_WRITE","SEQUENTIAL_READ")
472      Result='NO'
473   CASE ("OPEN_READ", "OPEN_WRITE", "OPEN_COMMIT_WRITE")
474      Result='REQUIRE'
475   CASE ("OPEN_COMMIT_READ","PARALLEL_IO")
476      Result='NO'
477   CASE ("SELF_DESCRIBING","SUPPORT_METADATA","SUPPORT_3D_FIELDS")
478      Result='YES'
479   CASE ("MEDIUM")
480      Result ='FILE'
481   CASE DEFAULT
482      Result = 'No Result for that inquiry!'
483   END SELECT
484   Status=WRF_NO_ERR
485   return
486 end subroutine ext_gr1_inquiry
488 !*****************************************************************************
490 SUBROUTINE ext_gr1_inquire_opened ( DataHandle, FileName , FileStat, Status )
492   USE gr1_data_info
493   IMPLICIT NONE
494 #include "wrf_status_codes.h"
495 #include "wrf_io_flags.h"
496   INTEGER ,       INTENT(IN)  :: DataHandle
497   CHARACTER*(*) :: FileName
498   INTEGER ,       INTENT(OUT) :: FileStat
499   INTEGER ,       INTENT(OUT) :: Status
501   call wrf_debug ( DEBUG , 'Entering ext_gr1_inquire_opened')
503   FileStat = WRF_NO_ERR
504   if ((DataHandle .ge. firstFileHandle) .and. &
505        (DataHandle .le. maxFileHandles)) then
506      FileStat = FileStatus(DataHandle)
507   else
508      FileStat = WRF_FILE_NOT_OPENED
509   endif
510   
511   Status = FileStat
513   RETURN
514 END SUBROUTINE ext_gr1_inquire_opened
516 !*****************************************************************************
518 SUBROUTINE ext_gr1_ioclose ( DataHandle, Status )
520   USE gr1_data_info
521   IMPLICIT NONE
522 #include "wrf_status_codes.h"
523   INTEGER DataHandle, Status
524   INTEGER istat
525   INTEGER ierr
526   character(len=1000) :: outstring
527   character :: lf
528   lf=char(10)
529      
530   call wrf_debug ( DEBUG , 'Entering ext_gr1_ioclose')
532   Status = WRF_NO_ERR
534   CALL write_file(FileFd(DataHandle), lf//'<METADATA>'//lf,ierr)
535   outstring = &
536        '<!-- The following are fields that were supplied to the WRF I/O API.'//lf//&
537        'Many variables (but not all) are redundant with the variables within '//lf//&
538        'the grib headers.  They are stored here, as METADATA, so that the '//lf//&
539        'WRF I/O API has simple access to these variables.-->'
540   CALL write_file(FileFd(DataHandle), trim(outstring), ierr)
541   if (trim(ti_output(DataHandle)) /= '') then
542      CALL write_file(FileFd(DataHandle), trim(ti_output(DataHandle)), ierr)
543      CALL write_file(FileFd(DataHandle), lf, ierr)
544   endif
545   if (trim(td_output(DataHandle)) /= '') then
546      CALL write_file(FileFd(DataHandle), trim(td_output(DataHandle)), ierr)
547      CALL write_file(FileFd(DataHandle), lf, ierr)
548   endif
549   CALL write_file(FileFd(DataHandle), '</METADATA>'//lf,ierr)
550   ti_output(DataHandle) = ''
551   td_output(DataHandle) = ''
552   if (ierr .ne. 0) then
553      Status = WRF_WARN_WRITE_RONLY_FILE
554   endif
555   CALL close_file(FileFd(DataHandle))
557   used(DataHandle) = .false.
559   RETURN
560 END SUBROUTINE ext_gr1_ioclose
562 !*****************************************************************************
564 SUBROUTINE ext_gr1_write_field( DataHandle , DateStrIn , VarName , &
565      Field , FieldType , Comm , IOComm, &
566      DomainDesc , MemoryOrder , Stagger , &
567      DimNames , &
568      DomainStart , DomainEnd , &
569      MemoryStart , MemoryEnd , &
570      PatchStart , PatchEnd , &
571      Status )
573   USE gr1_data_info
574   IMPLICIT NONE
575 #include "wrf_status_codes.h"
576 #include "wrf_io_flags.h"
577   INTEGER ,       INTENT(IN)    :: DataHandle 
578   CHARACTER*(*) :: DateStrIn
579   CHARACTER(DateStrLen) :: DateStr
580   CHARACTER*(*) :: VarName
581   CHARACTER*120 :: OutName
582   CHARACTER(120) :: TmpVarName
583   integer                       ,intent(in)    :: FieldType
584   integer                       ,intent(inout) :: Comm
585   integer                       ,intent(inout) :: IOComm
586   integer                       ,intent(in)    :: DomainDesc
587   character*(*)                 ,intent(in)    :: MemoryOrder
588   character*(*)                 ,intent(in)    :: Stagger
589   character*(*) , dimension (*) ,intent(in)    :: DimNames
590   integer ,dimension(*)         ,intent(in)    :: DomainStart, DomainEnd
591   integer ,dimension(*)         ,intent(in)    :: MemoryStart, MemoryEnd
592   integer ,dimension(*)         ,intent(in)    :: PatchStart,  PatchEnd
593   integer                       ,intent(out)   :: Status
594   integer                                      :: ierror
595   character (120)                         :: msg
596   integer :: xsize, ysize, zsize
597   integer :: x, y, z
598   integer :: x_start,x_end,y_start,y_end,z_start,z_end,ndim
599   integer :: idx
600   integer :: proj_center_flag
601   logical :: vert_stag = .false.
602   integer :: levelnum
603   real, DIMENSION(:,:), POINTER :: data,tmpdata
604   integer, DIMENSION(:), POINTER :: mold
605   integer :: istat
606   integer :: accum_period
607   integer :: size
608   integer, dimension(1000) :: level1, level2
609   real, DIMENSION( 1:1,MemoryStart(1):MemoryEnd(1), &
610                    MemoryStart(2):MemoryEnd(2), &
611                    MemoryStart(3):MemoryEnd(3) ) :: Field
612   real    :: fcst_secs
613   logical :: soil_layers, fraction
614   integer :: vert_unit
615   integer :: abc(2,2,2)
616   integer :: def(8)
617   logical :: output = .true.
618   integer :: idx1, idx2, idx3
619   integer :: this_domain
620   logical :: new_domain
621   real    :: region_center_lat, region_center_lon
622   integer :: dom_xsize, dom_ysize;
623   integer :: ierr
625   call wrf_debug ( DEBUG , 'Entering ext_gr1_write_field for parameter'//VarName)
627   !
628   ! If DateStr is all 0's, we reset it to StartDate (if StartDate exists).  
629   !   For some reason, 
630   !   in idealized simulations, StartDate is 0001-01-01_00:00:00 while
631   !   the first DateStr is 0000-00-00_00:00:00.  
632   !
633   if (DateStrIn .eq. '0000-00-00_00:00:00') then
634      if (StartDate .ne. '') then
635         DateStr = TRIM(StartDate)
636      else
637         DateStr = '0001-01-01_00:00:00'
638      endif
639   else
640      DateStr = DateStrIn
641   endif
643   !
644   ! Check if this is a domain that we haven't seen yet.  If so, add it to 
645   !   the list of domains.
646   !
647   this_domain = 0
648   new_domain = .false.
649   do idx = 1, max_domain
650      if (DomainDesc .eq. domains(idx)) then
651         this_domain = idx
652      endif
653   enddo
654   if (this_domain .eq. 0) then
655      max_domain = max_domain + 1
656      domains(max_domain) = DomainDesc
657      this_domain = max_domain
658      new_domain = .true.
659   endif
661   !
662   ! If the time has changed, we open a new file.  This is a kludge to get
663   !   around slowness in WRF that occurs when opening a new data file the
664   !   standard way.
665   !
666 #ifdef GRIB_ONE_TIME_PER_FILE
667   if (lastDateStr(this_domain) .ne. DateStr) then
668      write(DataFile(DataHandle),'(A8,i2.2,A1,A19)') 'wrfout_d',this_domain,'_',DateStr
669      print *,'Opening new file: ',trim(DataFile(DataHandle))
671      call ext_gr1_ioclose ( DataHandle, Status )
672      CALL open_file(TRIM(DataFile(DataHandle)), 'w', FileFd(DataHandle), ierr)
673      if (ierr .ne. 0) then
674         print *,'Could not open new file: ',DataFile(DataHandle)
675         print *,'  Appending to old file.'
676      else
677         ! Just set used back to .true. here, since ioclose set it to false.
678         used(DataHandle) = .true.
679      endif
680      td_output(DataHandle) = ''
681   endif
682   lastDateStr(this_domain) = DateStr
683 #endif
685   output = .true.
686   zsize = 1
687   xsize = 1
688   ysize = 1
689   OutName = VarName
690   soil_layers = .false.
691   fraction = .false.
693   ! First, handle then special cases for the boundary data.
695   CALL get_dims(MemoryOrder, PatchStart, PatchEnd, ndim, x_start, x_end, &
696        y_start, y_end,z_start,z_end)
697   xsize = x_end - x_start + 1
698   ysize = y_end - y_start + 1
699   zsize = z_end - z_start + 1
701   do idx = 1, len(MemoryOrder)
702      if ((MemoryOrder(idx:idx) .eq. 'Z') .and. &
703           (DimNames(idx) .eq. 'soil_layers_stag')) then
704         soil_layers = .true.
705      else if ((OutName .eq. 'LANDUSEF') .or. (OutName .eq. 'SOILCBOT') .or. &
706           (OutName .eq. 'SOILCTOP')) then
707         fraction = .true.
708      endif
709   enddo
711   if (.not. ASSOCIATED(grid_info)) then
712      CALL get_grid_info_size(size)
713      ALLOCATE(grid_info(1:size), STAT=istat)
714      if (istat .eq. -1) then
715         DEALLOCATE(grid_info)
716         Status = WRF_ERR_FATAL_BAD_FILE_STATUS
717         return
718      endif
719   endif
720      
722   if (new_domain) then
723      ALLOCATE(firstdata(this_domain)%snod(xsize,ysize))
724      firstdata(this_domain)%snod(:,:) = 0.0
725      ALLOCATE(lastdata(this_domain)%rainc(xsize,ysize))
726      lastdata(this_domain)%rainc(:,:) = 0.0
727      ALLOCATE(lastdata(this_domain)%rainnc(xsize,ysize))
728      lastdata(this_domain)%rainnc(:,:) = 0.0
729   endif
731   if (zsize .eq. 0) then 
732      zsize = 1
733   endif
735   ALLOCATE(data(1:xsize,1:ysize), STAT=istat)
736   ALLOCATE(mold(1:ysize), STAT=istat)
737   ALLOCATE(tmpdata(1:xsize,1:ysize), STAT=istat)
739   if (OutName .eq. 'ZNU') then
740      do idx = 1, zsize
741         half_eta(idx) = Field(1,idx,1,1)
742      enddo
743   endif
745   if (OutName .eq. 'ZNW') then
746      do idx = 1, zsize
747         full_eta(idx) = Field(1,idx,1,1)
748      enddo
749   endif
751   if (OutName .eq. 'ZS') then
752      do idx = 1, zsize
753         soil_depth(idx) = Field(1,idx,1,1)
754      enddo
755   endif
757   if (OutName .eq. 'DZS') then
758      do idx = 1, zsize
759         soil_thickness(idx) = Field(1,idx,1,1)
760      enddo
761   endif
764   if ((xsize .lt. 1) .or. (ysize .lt. 1)) then
765      write(msg,*) 'Cannot output field with memory order: ', &
766           MemoryOrder,Varname
767      call wrf_message(msg)
768      return
769   endif
770      
771   call get_vert_stag(OutName,Stagger,vert_stag)
773   do idx = 1, zsize
774      call gr1_get_levels(OutName, idx, zsize, soil_layers, vert_stag, fraction, &
775           vert_unit, level1(idx), level2(idx))
776   enddo
778   ! 
779   ! Get the center lat/lon for the area being output.  For some cases (such
780   !    as for boundary areas, the center of the area is different from the
781   !    center of the model grid.
782   !
783   if (index(Stagger,'X') .le. 0) then
784      dom_xsize = full_xsize - 1
785   else
786      dom_xsize = full_xsize
787   endif
788   if (index(Stagger,'Y') .le. 0) then
789      dom_ysize = full_ysize - 1
790   else
791      dom_ysize = full_ysize
792   endif
794   CALL get_region_center(MemoryOrder, projection, center_lat, center_lon, &
795        dom_xsize, dom_ysize, dx, dy, proj_central_lon, proj_center_flag, &
796        truelat1, truelat2, xsize, ysize, region_center_lat, region_center_lon)
798   if ( .not. opened(DataHandle)) then
799      Status = WRF_WARN_FILE_NOT_OPENED
800      return
801   endif
804   if (opened(DataHandle) .and. committed(DataHandle)) then
807 #ifdef OUTPUT_FULL_PRESSURE
809      ! 
810      ! The following is a kludge to output full pressure instead of the two 
811      !  fields of base-state pressure and pressure perturbation.
812      !
813      ! This code can be turned on by adding -DOUTPUT_FULL_PRESSURE to the 
814      !  compile line
815      !
816      
817      if ((OutName .eq. 'P') .or. (OutName.eq.'PB')) then
818         do idx = 1, len(MemoryOrder)
819             if (MemoryOrder(idx:idx) .eq. 'X') then
820               idx1=idx
821            endif
822            if (MemoryOrder(idx:idx) .eq. 'Y') then
823               idx2=idx
824            endif
825            if (MemoryOrder(idx:idx) .eq. 'Z') then
826               idx3=idx
827            endif
828         enddo
830         ! 
831         ! Allocate space for pressure values (this variable holds 
832         !   base-state pressure or pressure perturbation to be used 
833         !   later to sum base-state and perturbation pressure to get full 
834         !   pressure).
835         !
837         if (.not. ASSOCIATED(pressure(this_domain)%vals)) then
838            ALLOCATE(pressure(this_domain)%vals(MemoryStart(1):MemoryEnd(1), &
839                 MemoryStart(2):MemoryEnd(2),MemoryStart(3):MemoryEnd(3)))
840         endif
841         if (DateStr .NE. &
842              pressure(this_domain)%lastDateStr) then
843            pressure(this_domain)%newtime = .true.
844         endif
845         if (pressure(this_domain)%newtime) then
846            pressure(this_domain)%vals = Field(1,:,:,:)
847            pressure(this_domain)%newtime = .false.
848            output = .false.
849         else 
850            output = .true.
851         endif
852         pressure(this_domain)%lastDateStr=DateStr
853      endif
854 #endif
856 #ifdef OUTPUT_FULL_GEOPOTENTIAL
858      ! 
859      ! The following is a kludge to output full geopotential height instead 
860      !  of the two fields of base-state geopotential and perturbation 
861      !  geopotential.
862      !
863      ! This code can be turned on by adding -DOUTPUT_FULL_GEOPOTENTIAL to the 
864      !  compile line
865      !
866      
867      if ((OutName .eq. 'PHB') .or. (OutName.eq.'PH')) then
868         do idx = 1, len(MemoryOrder)
869             if (MemoryOrder(idx:idx) .eq. 'X') then
870               idx1=idx
871            endif
872            if (MemoryOrder(idx:idx) .eq. 'Y') then
873               idx2=idx
874            endif
875            if (MemoryOrder(idx:idx) .eq. 'Z') then
876               idx3=idx
877            endif
878         enddo
880         ! 
881         ! Allocate space for geopotential values (this variable holds 
882         !   geopotential to be used 
883         !   later to sum base-state and perturbation to get full 
884         !   geopotential).
885         !
887         if (.not. ASSOCIATED(geopotential(this_domain)%vals)) then
888            ALLOCATE(geopotential(this_domain)%vals(MemoryStart(1):MemoryEnd(1), &
889                 MemoryStart(2):MemoryEnd(2),MemoryStart(3):MemoryEnd(3)))
890         endif
891         if (DateStr .NE. &
892              geopotential(this_domain)%lastDateStr) then
893            geopotential(this_domain)%newtime = .true.
894         endif
895         if (geopotential(this_domain)%newtime) then
896            geopotential(this_domain)%vals = Field(1,:,:,:)
897            geopotential(this_domain)%newtime = .false.
898            output = .false.
899         else 
900            output = .true.
901         endif
902         geopotential(this_domain)%lastDateStr=DateStr
903      endif
904 #endif
906      if (output) then 
907         if (StartDate == '') then
908            StartDate = DateStr
909         endif
910         CALL geth_idts(DateStr,StartDate,fcst_secs)
911         
912         if (center_lat .lt. 0) then
913            proj_center_flag = 2
914         else
915            proj_center_flag = 1
916         endif
917          
918         do z = 1, zsize
919            SELECT CASE (MemoryOrder)
920            CASE ('XYZ')
921               data = Field(1,1:xsize,1:ysize,z)
922            CASE ('XZY')
923               data = Field(1,1:xsize,z,1:ysize)
924            CASE ('YXZ')
925               do x = 1,xsize
926                  do y = 1,ysize
927                     data(x,y) = Field(1,y,x,z)
928                  enddo
929               enddo
930            CASE ('YZX')
931               do x = 1,xsize
932                  do y = 1,ysize
933                     data(x,y) = Field(1,y,z,x)
934                  enddo
935               enddo
936            CASE ('ZXY')
937               data = Field(1,z,1:xsize,1:ysize)
938            CASE ('ZYX')
939               do x = 1,xsize
940                  do y = 1,ysize
941                     data(x,y) = Field(1,z,y,x)
942                  enddo
943               enddo
944            CASE ('XY')
945               data = Field(1,1:xsize,1:ysize,1)
946            CASE ('YX')
947               do x = 1,xsize
948                  do y = 1,ysize
949                     data(x,y) = Field(1,y,x,1)
950                  enddo
951               enddo
953            CASE ('XSZ')
954               do x = 1,xsize
955                  do y = 1,ysize
956                     data(x,y) = Field(1,y,z,x)
957                  enddo
958               enddo
959            CASE ('XEZ')
960               do x = 1,xsize
961                  do y = 1,ysize
962                     data(x,y) = Field(1,y,z,x)
963                  enddo
964               enddo
965            CASE ('YSZ')
966               do x = 1,xsize
967                  do y = 1,ysize
968                     data(x,y) = Field(1,x,z,y)
969                  enddo
970               enddo
971            CASE ('YEZ')
972               do x = 1,xsize
973                  do y = 1,ysize
974                     data(x,y) = Field(1,x,z,y)
975                  enddo
976               enddo
978            CASE ('XS')
979               do x = 1,xsize
980                  do y = 1,ysize
981                     data(x,y) = Field(1,y,x,1)
982                  enddo
983               enddo
984            CASE ('XE')
985               do x = 1,xsize
986                  do y = 1,ysize
987                     data(x,y) = Field(1,y,x,1)
988                  enddo
989               enddo
990            CASE ('YS')
991               do x = 1,xsize
992                  do y = 1,ysize
993                     data(x,y) = Field(1,x,y,1)
994                  enddo
995               enddo
996            CASE ('YE')
997               do x = 1,xsize
998                  do y = 1,ysize
999                     data(x,y) = Field(1,x,y,1)
1000                  enddo
1001               enddo
1003            CASE ('Z')
1004               data(1,1) = Field(1,z,1,1)
1005            CASE ('z')
1006               data(1,1) = Field(1,z,1,1)
1007            CASE ('C')
1008               data = Field(1,1:xsize,1:ysize,z)
1009            CASE ('c')
1010               data = Field(1,1:xsize,1:ysize,z)
1011            CASE ('0')
1012               data(1,1) = Field(1,1,1,1)
1013            END SELECT
1015            ! 
1016            ! Here, we convert any integer fields to real
1017            !
1018            if (FieldType == WRF_INTEGER) then
1019               mold = 0
1020               do idx=1,xsize
1021                  !
1022                  ! The parentheses around data(idx,:) are needed in order
1023                  !   to fix a bug with transfer with the xlf compiler on NCAR's
1024                  !   IBM (bluesky).
1025                  !
1026                  data(idx,:)=transfer((data(idx,:)),mold)
1027               enddo
1028            endif
1029            ! 
1030            ! Here, we do any necessary conversions to the data.
1031            !
1032            
1033            ! Potential temperature is sometimes passed in as perturbation 
1034            !   potential temperature (i.e., POT-300).  Other times (i.e., from 
1035            !   WRF SI), it is passed in as full potential temperature.
1036            ! Here, we convert to full potential temperature by adding 300
1037            !   only if POT < 200 K.
1038            !
1039            if (OutName == 'T') then
1040               if (data(1,1) < 200) then
1041                  data = data + 300
1042               endif
1043            endif
1045            ! 
1046            ! For precip, we setup the accumulation period, and output a precip
1047            !    rate for time-step precip.
1048            !
1049            if (OutName .eq. 'RAINNCV') then
1050               ! Convert time-step precip to precip rate.
1051               data = data/timestep
1052               accum_period = 0
1053            else
1054               accum_period = 0
1055            endif
1057 #ifdef OUTPUT_FULL_PRESSURE
1058            !
1059            ! Computation of full-pressure off by default since there are 
1060            !  uses for base-state and perturbation (i.e., restarts
1061            !
1062             if ((OutName .eq. 'P') .or. (OutName.eq.'PB')) then
1063                if (idx3 .eq. 1) then
1064                   data = data + &
1065                        pressure(this_domain)%vals(z, &
1066                        patchstart(2):patchend(2),patchstart(3):patchend(3))
1067                elseif (idx3 .eq. 2) then
1068                   data = data + &
1069                        pressure(this_domain)%vals(patchstart(1):patchend(1), &
1070                        z,patchstart(3):patchend(3))
1071                elseif (idx3 .eq. 3) then
1072                   data = data + &
1073                        pressure(this_domain)%vals(patchstart(1):patchend(1), &
1074                        patchstart(2):patchend(2),z)
1075                else
1076                   call wrf_message ('error in idx3, continuing')
1077                endif
1079                OutName = 'P'
1080             endif
1081 #endif
1083 #ifdef OUTPUT_FULL_GEOPOTENTIAL
1084            !
1085            ! Computation of full-geopotential off by default since there are 
1086            !  uses for base-state and perturbation (i.e., restarts
1087            !
1088             if ((OutName .eq. 'PHB') .or. (OutName.eq.'PH')) then
1089                if (idx3 .eq. 1) then
1090                   data = data + &
1091                        geopotential(this_domain)%vals(z, &
1092                        patchstart(2):patchend(2),patchstart(3):patchend(3))
1093                elseif (idx3 .eq. 2) then
1094                   data = data + &
1095                        geopotential(this_domain)%vals(patchstart(1):patchend(1), &
1096                        z,patchstart(3):patchend(3))
1097                elseif (idx3 .eq. 3) then
1098                   data = data + &
1099                        geopotential(this_domain)%vals(patchstart(1):patchend(1), &
1100                        patchstart(2):patchend(2),z)
1101                else
1102                   call wrf_message ('error in idx3, continuing')
1103                endif
1105                OutName = 'PHP'
1106             endif
1107 #endif
1109            !
1110            !    Output current level
1111            !
1112            CALL load_grid_info(OutName, StartDate, vert_unit, level1(z), &
1113                 level2(z), fcst_secs, accum_period, wg_grid_id, projection, &
1114                 xsize, ysize, region_center_lat, region_center_lon, dx, dy, &
1115                 proj_central_lon, proj_center_flag, truelat1, truelat2, &
1116                 grib_tables, grid_info)
1117            
1118            !
1119            ! Here, we copy data to a temporary array.  After write_grib,
1120            !    we copy back from the temporary array to the permanent
1121            !    array.  write_grib modifies data.  For certain fields that
1122            !    we use below, we want the original (unmodified) data 
1123            !    values.  This kludge assures that we have the original
1124            !    values.
1125            !
1127            if ((OutName .eq. 'RAINC') .or. (OutName .eq. 'RAINNC') .or. &
1128                 (OutName .eq. 'SNOWH')) then
1129               tmpdata(:,:) = data(:,:)
1130            endif
1132            CALL write_grib(grid_info, FileFd(DataHandle), data)
1134            if ((OutName .eq. 'RAINC') .or. (OutName .eq. 'RAINNC') .or. &
1135                 (OutName .eq. 'SNOWH')) then
1136               data(:,:) = tmpdata(:,:)
1137            endif
1139            CALL free_grid_info(grid_info)
1140            
1141            !
1142            ! If this is the total accumulated rain, call write_grib again 
1143            !   to output the accumulation since the last output time as well.
1144            !   This is somewhat of a kludge to meet the requirements of PF.
1145            !
1146            if ((OutName .eq. 'RAINC') .or. (OutName .eq. 'RAINNC') .or. &
1147                 (OutName .eq. 'SNOWH')) then
1148               if (OutName .eq. 'RAINC') then
1149                  data(:,:) = data(:,:) - lastdata(this_domain)%rainc(:,:)
1150                  lastdata(this_domain)%rainc(:,:) = tmpdata(:,:)
1151                  accum_period = fcst_secs - &
1152                       lastdata(this_domain)%fcst_secs_rainc
1153                  lastdata(this_domain)%fcst_secs_rainc = fcst_secs
1154                  TmpVarName = 'ACPCP'
1155               else if (OutName .eq. 'RAINNC') then
1156                  tmpdata(:,:) = data(:,:)
1157                  data(:,:) = data(:,:) - lastdata(this_domain)%rainnc(:,:)
1158                  lastdata(this_domain)%rainnc(:,:) = tmpdata(:,:)
1159                  accum_period = fcst_secs - &
1160                       lastdata(this_domain)%fcst_secs_rainnc
1161                  lastdata(this_domain)%fcst_secs_rainnc = fcst_secs
1162                  TmpVarName = 'NCPCP'
1163               else if (OutName .eq. 'SNOWH') then
1164                  if (fcst_secs .eq. 0) then
1165                     firstdata(this_domain)%snod(:,:) = data(:,:)
1166                  endif
1167                  data(:,:) = data(:,:) - firstdata(this_domain)%snod(:,:)
1168                  TmpVarName = 'SNOWCU'
1169               endif
1171               CALL load_grid_info(TmpVarName, StartDate, vert_unit, level1(z),&
1172                    level2(z), fcst_secs, accum_period, wg_grid_id, &
1173                    projection, xsize, ysize, region_center_lat, &
1174                    region_center_lon, dx, dy, proj_central_lon, &
1175                    proj_center_flag, truelat1, truelat2, grib_tables, &
1176                    grid_info)
1177            
1178               CALL write_grib(grid_info, FileFd(DataHandle), data)
1179               CALL free_grid_info(grid_info)
1180            endif
1182         enddo
1183      endif
1184   endif
1186   deallocate(data, STAT = istat)
1187   deallocate(mold, STAT = istat)
1188   deallocate(tmpdata, STAT = istat)
1190   Status = WRF_NO_ERR
1192   call wrf_debug ( DEBUG , 'Leaving ext_gr1_write_field')
1194   RETURN
1195 END SUBROUTINE ext_gr1_write_field
1197 !*****************************************************************************
1199 SUBROUTINE ext_gr1_read_field ( DataHandle , DateStr , VarName , Field , &
1200      FieldType , Comm , IOComm, DomainDesc , MemoryOrder , Stagger ,     &
1201      DimNames , DomainStart , DomainEnd , MemoryStart , MemoryEnd ,      &
1202      PatchStart , PatchEnd ,  Status )
1204   USE gr1_data_info
1205   IMPLICIT NONE  
1206 #include "wrf_status_codes.h"
1207 #include "wrf_io_flags.h"
1208   INTEGER ,       INTENT(IN)    :: DataHandle 
1209   CHARACTER*(*) :: DateStr
1210   CHARACTER*(*) :: VarName
1211   CHARACTER (len=400) :: msg
1212   integer                       ,intent(inout)    :: FieldType
1213   integer                       ,intent(inout)    :: Comm
1214   integer                       ,intent(inout)    :: IOComm
1215   integer                       ,intent(inout)    :: DomainDesc
1216   character*(*)                 ,intent(inout)    :: MemoryOrder
1217   character*(*)                 ,intent(inout)    :: Stagger
1218   character*(*) , dimension (*) ,intent(inout)    :: DimNames
1219   integer ,dimension(*)         ,intent(inout)    :: DomainStart, DomainEnd
1220   integer ,dimension(*)         ,intent(inout)    :: MemoryStart, MemoryEnd
1221   integer ,dimension(*)         ,intent(inout)    :: PatchStart,  PatchEnd
1222   integer                       ,intent(out)      :: Status
1223   INTEGER                       ,intent(out)      :: Field(*)
1224   integer   :: ndim,x_start,x_end,y_start,y_end,z_start,z_end
1225   integer   :: zidx
1226   REAL, DIMENSION(:,:), POINTER :: data
1227   logical                     :: vert_stag
1228   logical                     :: soil_layers
1229   integer                     :: level1,level2
1231   integer                     :: parmid
1232   integer                     :: vert_unit
1233   integer                     :: grb_index
1234   integer                     :: numcols, numrows
1235   integer                     :: data_allocated
1236   integer                     :: istat
1237   integer                     :: tablenum
1238   integer                     :: di
1239   integer                     :: last_grb_index
1241   call wrf_debug ( DEBUG , 'Entering ext_gr1_read_field')
1243   !
1244   ! Get dimensions of data.  
1245   ! Assume that the domain size in the input data is the same as the Domain 
1246   !     Size from the input arguments.
1247   !
1248   
1249   CALL get_dims(MemoryOrder,DomainStart,DomainEnd,ndim,x_start,x_end,y_start, &
1250        y_end,z_start,z_end) 
1252   !
1253   ! Get grib parameter id
1254   !
1255   CALL GET_GRIB_PARAM(grib_tables, VarName, center, subcenter, parmtbl, &
1256        tablenum, parmid)
1258   !
1259   ! Setup the vertical unit and levels
1260   !
1261   CALL get_vert_stag(VarName,Stagger,vert_stag)
1262   CALL get_soil_layers(VarName,soil_layers)
1264   !
1265   ! Loop over levels, grabbing data from each level, then assembling into a 
1266   !   3D array.
1267   !
1268   data_allocated = 0
1269   last_grb_index = -1
1270   do zidx = z_start,z_end
1271      
1272      CALL gr1_get_levels(VarName,zidx,z_end-z_start,soil_layers,vert_stag, &
1273           .false., vert_unit,level1,level2)
1274      
1275      CALL GET_GRIB_INDEX_VALIDTIME_GUESS(fileinfo(DataHandle)%fileindex(:), center, &
1276           subcenter, parmtbl, parmid,DateStr,vert_unit,level1, &
1277           level2, last_grb_index + 1, grb_index)
1278      if (grb_index < 0) then
1279         write(msg,*)'Field not found: parmid: ',VarName,parmid,DateStr, &
1280              vert_unit,level1,level2
1281         call wrf_debug (DEBUG , msg)
1282         cycle
1283      endif
1285      if (data_allocated .eq. 0) then
1286         CALL GET_SIZEOF_GRID(fileinfo(DataHandle)%fileindex(:),grb_index,numcols,numrows)
1287         allocate(data(z_start:z_end,1:numcols*numrows),stat=istat)
1288         data_allocated = 1
1289      endif
1291      CALL READ_GRIB(fileinfo(DataHandle)%fileindex(:), FileFd(DataHandle), grb_index, &
1292           data(zidx,:))
1294      !
1295      ! Transpose data into the order specified by MemoryOrder, setting only 
1296      !   entries within the memory dimensions
1297      !
1298      CALL get_dims(MemoryOrder, MemoryStart, MemoryEnd, ndim, x_start, x_end, &
1299           y_start, y_end,z_start,z_end)
1301      if(FieldType == WRF_DOUBLE) then
1302         di = 2
1303      else 
1304         di = 1
1305      endif
1307      ! 
1308      ! Here, we do any necessary conversions to the data.
1309      !
1310      ! The WRF executable (wrf.exe) expects perturbation potential
1311      !   temperature.  However, real.exe expects full potential T.
1312      ! So, if the program is WRF, subtract 300 from Potential Temperature 
1313      !   to get perturbation potential temperature.
1314      !
1315      if (VarName == 'T') then
1316         if ( &
1317              (InputProgramName .eq. 'REAL_EM') .or. &
1318              (InputProgramName .eq. 'IDEAL') .or. &
1319              (InputProgramName .eq. 'NDOWN_EM')) then
1320            data(zidx,:) = data(zidx,:) - 300
1321         endif
1322      endif
1324      CALL Transpose_grib(MemoryOrder, di, FieldType, Field, &
1325           MemoryStart(1), MemoryEnd(1), MemoryStart(2), MemoryEnd(2), &
1326           MemoryStart(3), MemoryEnd(3), &
1327           data(zidx,:), zidx, numrows, numcols)
1329      if (zidx .eq. z_end) then
1330         data_allocated = 0
1331         deallocate(data)
1332      endif
1334      last_grb_index = grb_index
1336   enddo
1338   Status = WRF_NO_ERR
1339   if (grb_index < 0) Status = WRF_WARN_VAR_NF
1340   call wrf_debug ( DEBUG , 'Leaving ext_gr1_read_field')
1342   RETURN
1343 END SUBROUTINE ext_gr1_read_field
1345 !*****************************************************************************
1347 SUBROUTINE ext_gr1_get_next_var ( DataHandle, VarName, Status )
1349   USE gr1_data_info
1350   IMPLICIT NONE
1351 #include "wrf_status_codes.h"
1352   INTEGER ,       INTENT(IN)  :: DataHandle
1353   CHARACTER*(*) :: VarName
1354   INTEGER ,       INTENT(OUT) :: Status
1356   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_next_var')
1358   call wrf_message ( 'WARNING: ext_gr1_get_next_var is not supported.')
1360   Status = WRF_WARN_NOOP
1362   RETURN
1363 END SUBROUTINE ext_gr1_get_next_var
1365 !*****************************************************************************
1367 subroutine ext_gr1_end_of_frame(DataHandle, Status)
1369   USE gr1_data_info
1370   implicit none
1371 #include "wrf_status_codes.h"
1372   integer               ,intent(in)     :: DataHandle
1373   integer               ,intent(out)    :: Status
1375   call wrf_debug ( DEBUG , 'Entering ext_gr1_end_of_frame')
1377   Status = WRF_WARN_NOOP
1379   return
1380 end subroutine ext_gr1_end_of_frame
1382 !*****************************************************************************
1384 SUBROUTINE ext_gr1_iosync ( DataHandle, Status )
1386   USE gr1_data_info  
1387   IMPLICIT NONE
1388 #include "wrf_status_codes.h"
1389   INTEGER ,       INTENT(IN)  :: DataHandle
1390   INTEGER ,       INTENT(OUT) :: Status
1392   call wrf_debug ( DEBUG , 'Entering ext_gr1_iosync')
1394   Status = WRF_NO_ERR
1395   if (DataHandle .GT. 0) then
1396      CALL flush_file(FileFd(DataHandle))
1397   else
1398      Status = WRF_WARN_TOO_MANY_FILES
1399   endif
1401   RETURN
1402 END SUBROUTINE ext_gr1_iosync
1404 !*****************************************************************************
1406 SUBROUTINE ext_gr1_inquire_filename ( DataHandle, FileName , FileStat, &
1407      Status )
1409   USE gr1_data_info
1410   IMPLICIT NONE
1411 #include "wrf_status_codes.h"
1412 #include "wrf_io_flags.h"
1413   INTEGER ,       INTENT(IN)  :: DataHandle
1414   CHARACTER*(*) :: FileName
1415   INTEGER ,       INTENT(OUT) :: FileStat
1416   INTEGER ,       INTENT(OUT) :: Status
1417   CHARACTER *80   SysDepInfo
1419   call wrf_debug ( DEBUG , 'Entering ext_gr1_inquire_filename')
1421   FileName = DataFile(DataHandle) 
1423   if ((DataHandle .ge. firstFileHandle) .and. &
1424        (DataHandle .le. maxFileHandles)) then
1425      FileStat = FileStatus(DataHandle)
1426   else
1427      FileStat = WRF_FILE_NOT_OPENED
1428   endif
1429   
1430   Status = WRF_NO_ERR
1432   RETURN
1433 END SUBROUTINE ext_gr1_inquire_filename
1435 !*****************************************************************************
1437 SUBROUTINE ext_gr1_get_var_info ( DataHandle , VarName , NDim , &
1438      MemoryOrder , Stagger , DomainStart , DomainEnd , WrfType, Status )
1440   USE gr1_data_info
1441   IMPLICIT NONE
1442 #include "wrf_status_codes.h"
1443   integer               ,intent(in)     :: DataHandle
1444   character*(*)         ,intent(in)     :: VarName
1445   integer               ,intent(out)    :: NDim
1446   character*(*)         ,intent(out)    :: MemoryOrder
1447   character*(*)         ,intent(out)    :: Stagger
1448   integer ,dimension(*) ,intent(out)    :: DomainStart, DomainEnd
1449   integer               ,intent(out)    :: WrfType
1450   integer               ,intent(out)    :: Status
1452   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_info')
1454   CALL wrf_message('ext_gr1_get_var_info not supported for grib version1 data')
1455   Status = WRF_NO_ERR
1457   RETURN
1458 END SUBROUTINE ext_gr1_get_var_info
1460 !*****************************************************************************
1462 SUBROUTINE ext_gr1_set_time ( DataHandle, DateStr, Status )
1464   USE gr1_data_info
1465   IMPLICIT NONE
1466 #include "wrf_status_codes.h"
1467   INTEGER ,       INTENT(IN)  :: DataHandle
1468   CHARACTER*(*) :: DateStr
1469   INTEGER ,       INTENT(OUT) :: Status
1470   integer       :: found_time
1471   integer       :: idx
1473   call wrf_debug ( DEBUG , 'Entering ext_gr1_set_time')
1475   found_time = 0
1476   do idx = 1,fileinfo(DataHandle)%NumberTimes
1477      if (fileinfo(DataHandle)%Times(idx) == DateStr) then
1478         found_time = 1
1479         fileinfo(DataHandle)%CurrentTime = idx
1480      endif
1481   enddo
1482   if (found_time == 0) then 
1483      Status = WRF_WARN_TIME_NF
1484   else
1485      Status = WRF_NO_ERR
1486   endif
1488   RETURN
1489 END SUBROUTINE ext_gr1_set_time
1491 !*****************************************************************************
1493 SUBROUTINE ext_gr1_get_next_time ( DataHandle, DateStr, Status )
1495   USE gr1_data_info
1496   IMPLICIT NONE
1497 #include "wrf_status_codes.h"
1498   INTEGER ,       INTENT(IN)  :: DataHandle
1499   CHARACTER*(*) , INTENT(OUT) :: DateStr
1500   INTEGER ,       INTENT(OUT) :: Status
1502   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_next_time')
1504   if (fileinfo(DataHandle)%CurrentTime == fileinfo(DataHandle)%NumberTimes) then
1505      Status = WRF_WARN_TIME_EOF
1506   else
1507      fileinfo(DataHandle)%CurrentTime = fileinfo(DataHandle)%CurrentTime + 1
1508      DateStr = fileinfo(DataHandle)%Times(fileinfo(DataHandle)%CurrentTime)
1509      Status = WRF_NO_ERR
1510   endif
1512   RETURN
1513 END SUBROUTINE ext_gr1_get_next_time
1515 !*****************************************************************************
1517 SUBROUTINE ext_gr1_get_previous_time ( DataHandle, DateStr, Status )
1519   USE gr1_data_info
1520   IMPLICIT NONE
1521 #include "wrf_status_codes.h"
1522   INTEGER ,       INTENT(IN)  :: DataHandle
1523   CHARACTER*(*) :: DateStr
1524   INTEGER ,       INTENT(OUT) :: Status
1526   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_previous_time')
1528   if (fileinfo(DataHandle)%CurrentTime <= 0) then
1529      Status = WRF_WARN_TIME_EOF
1530   else
1531      fileinfo(DataHandle)%CurrentTime = fileinfo(DataHandle)%CurrentTime - 1
1532      DateStr = fileinfo(DataHandle)%Times(fileinfo(DataHandle)%CurrentTime)
1533      Status = WRF_NO_ERR
1534   endif
1536   RETURN
1537 END SUBROUTINE ext_gr1_get_previous_time
1539 !******************************************************************************
1540 !* Start of get_var_ti_* routines
1541 !******************************************************************************
1543 SUBROUTINE ext_gr1_get_var_ti_real ( DataHandle,Element,  Varname, Data, &
1544      Count, Outcount, Status )
1546   USE gr1_data_info
1547   IMPLICIT NONE
1548 #include "wrf_status_codes.h"
1549   INTEGER ,       INTENT(IN)    :: DataHandle
1550   CHARACTER*(*) :: Element
1551   CHARACTER*(*) :: VarName 
1552   real ,          INTENT(OUT)   :: Data(*)
1553   INTEGER ,       INTENT(IN)    :: Count
1554   INTEGER ,       INTENT(OUT)   :: OutCount
1555   INTEGER ,       INTENT(OUT)   :: Status
1556   INTEGER          :: idx
1557   INTEGER          :: stat
1558   CHARACTER*(1000) :: VALUE
1560   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_ti_real')
1562   Status = WRF_NO_ERR
1563   
1564   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), "none", &
1565        Varname, Value, stat)
1566   if (stat /= 0) then
1567      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
1568      Status = WRF_WARN_VAR_NF
1569      RETURN
1570   endif
1572   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
1573   if (stat .ne. 0) then
1574      CALL wrf_message("Reading data from"//Value//"failed")
1575      Status = WRF_WARN_COUNT_TOO_LONG
1576      RETURN
1577   endif
1578   Outcount = idx
1580   RETURN
1581 END SUBROUTINE ext_gr1_get_var_ti_real 
1583 !*****************************************************************************
1585 SUBROUTINE ext_gr1_get_var_ti_real8 ( DataHandle,Element,  Varname, Data, &
1586      Count, Outcount, Status )
1588   USE gr1_data_info
1589   IMPLICIT NONE
1590 #include "wrf_status_codes.h"
1591   INTEGER ,       INTENT(IN)      :: DataHandle
1592   CHARACTER*(*) :: Element
1593   CHARACTER*(*) :: VarName 
1594   real*8 ,        INTENT(OUT)     :: Data(*)
1595   INTEGER ,       INTENT(IN)      :: Count
1596   INTEGER ,       INTENT(OUT)     :: OutCount
1597   INTEGER ,       INTENT(OUT)     :: Status
1598   INTEGER          :: idx
1599   INTEGER          :: stat
1600   CHARACTER*(1000) :: VALUE
1602   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_ti_real8')
1604   Status = WRF_NO_ERR
1605   
1606   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:),TRIM(Element),&
1607        "none",Varname,Value,stat)
1608   if (stat /= 0) then
1609      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
1610      Status = WRF_WARN_VAR_NF
1611      RETURN
1612   endif
1614   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
1615   if (stat .ne. 0) then
1616      CALL wrf_message("Reading data from"//Value//"failed")
1617      Status = WRF_WARN_COUNT_TOO_LONG
1618      RETURN
1619   endif
1620   Outcount = idx
1622   RETURN
1623 END SUBROUTINE ext_gr1_get_var_ti_real8 
1625 !*****************************************************************************
1627 SUBROUTINE ext_gr1_get_var_ti_double ( DataHandle,Element,  Varname, Data, &
1628      Count, Outcount, Status )
1629   USE gr1_data_info
1630   IMPLICIT NONE
1631 #include "wrf_status_codes.h"
1632   INTEGER ,       INTENT(IN)  :: DataHandle
1633   CHARACTER*(*) , INTENT(IN)  :: Element
1634   CHARACTER*(*) , INTENT(IN)  :: VarName
1635   real*8 ,            INTENT(OUT) :: Data(*)
1636   INTEGER ,       INTENT(IN)  :: Count
1637   INTEGER ,       INTENT(OUT)  :: OutCount
1638   INTEGER ,       INTENT(OUT) :: Status
1639   INTEGER          :: idx
1640   INTEGER          :: stat
1641   CHARACTER*(1000) :: VALUE
1643   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_ti_double')
1645   Status = WRF_NO_ERR
1646   
1647   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), &
1648        "none", Varname, &
1649        Value,stat)
1650   if (stat /= 0) then
1651      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
1652      Status = WRF_WARN_VAR_NF
1653      RETURN
1654   endif
1656   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
1657   if (stat .ne. 0) then
1658      CALL wrf_message("Reading data from"//Value//"failed")
1659      Status = WRF_WARN_COUNT_TOO_LONG
1660      RETURN
1661   endif
1662   Outcount = idx
1664   RETURN
1665 END SUBROUTINE ext_gr1_get_var_ti_double
1667 !*****************************************************************************
1669 SUBROUTINE ext_gr1_get_var_ti_integer ( DataHandle,Element,  Varname, Data, &
1670      Count, Outcount, Status )
1672   USE gr1_data_info
1673   IMPLICIT NONE
1674 #include "wrf_status_codes.h"
1675   INTEGER ,       INTENT(IN)       :: DataHandle
1676   CHARACTER*(*) :: Element
1677   CHARACTER*(*) :: VarName 
1678   integer ,       INTENT(OUT)      :: Data(*)
1679   INTEGER ,       INTENT(IN)       :: Count
1680   INTEGER ,       INTENT(OUT)      :: OutCount
1681   INTEGER ,       INTENT(OUT)      :: Status
1682   INTEGER          :: idx
1683   INTEGER          :: stat
1684   CHARACTER*(1000) :: VALUE
1686   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_ti_integer')
1688   Status = WRF_NO_ERR
1689   
1690   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), &
1691        "none", Varname, Value, stat)
1692   if (stat /= 0) then
1693      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
1694      Status = WRF_WARN_VAR_NF
1695      RETURN
1696   endif
1698   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
1699   if (stat .ne. 0) then
1700      CALL wrf_message("Reading data from"//Value//"failed")
1701      Status = WRF_WARN_COUNT_TOO_LONG
1702      RETURN
1703   endif
1704   Outcount = idx
1706   RETURN
1707 END SUBROUTINE ext_gr1_get_var_ti_integer 
1709 !*****************************************************************************
1711 SUBROUTINE ext_gr1_get_var_ti_logical ( DataHandle,Element,  Varname, Data, &
1712      Count, Outcount, Status )
1714   USE gr1_data_info
1715   IMPLICIT NONE
1716 #include "wrf_status_codes.h"
1717   INTEGER ,       INTENT(IN)       :: DataHandle
1718   CHARACTER*(*) :: Element
1719   CHARACTER*(*) :: VarName 
1720   logical ,       INTENT(OUT)      :: Data(*)
1721   INTEGER ,       INTENT(IN)       :: Count
1722   INTEGER ,       INTENT(OUT)      :: OutCount
1723   INTEGER ,       INTENT(OUT)      :: Status
1724   INTEGER          :: idx
1725   INTEGER          :: stat
1726   CHARACTER*(1000) :: VALUE
1728   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_ti_logical')
1730   Status = WRF_NO_ERR
1731   
1732   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), &
1733        "none", Varname, Value,stat)
1734   if (stat /= 0) then
1735      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
1736      Status = WRF_WARN_VAR_NF
1737      RETURN
1738   endif
1740   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
1741   if (stat .ne. 0) then
1742      CALL wrf_message("Reading data from"//Value//"failed")
1743      Status = WRF_WARN_COUNT_TOO_LONG
1744      RETURN
1745   endif
1746   Outcount = idx
1748   RETURN
1749 END SUBROUTINE ext_gr1_get_var_ti_logical 
1751 !*****************************************************************************
1753 SUBROUTINE ext_gr1_get_var_ti_char ( DataHandle,Element,  Varname, Data,  &
1754      Status )
1756   USE gr1_data_info
1757   IMPLICIT NONE
1758 #include "wrf_status_codes.h"
1759   INTEGER ,       INTENT(IN)  :: DataHandle
1760   CHARACTER*(*) :: Element
1761   CHARACTER*(*) :: VarName 
1762   CHARACTER*(*) :: Data
1763   INTEGER ,       INTENT(OUT) :: Status
1764   INTEGER       :: stat
1766   Status = WRF_NO_ERR
1767   
1768   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_ti_char')
1770   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), &
1771        "none", Varname, Data,stat)
1772   if (stat /= 0) then
1773      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
1774      Status = WRF_WARN_VAR_NF
1775      RETURN
1776   endif
1778   RETURN
1779 END SUBROUTINE ext_gr1_get_var_ti_char 
1781 !******************************************************************************
1782 !* End of get_var_ti_* routines
1783 !******************************************************************************
1786 !******************************************************************************
1787 !* Start of put_var_ti_* routines
1788 !******************************************************************************
1790 SUBROUTINE ext_gr1_put_var_ti_real ( DataHandle,Element,  Varname, Data, &
1791      Count,  Status )
1793   USE gr1_data_info
1794   IMPLICIT NONE
1795 #include "wrf_status_codes.h"
1796   INTEGER ,       INTENT(IN)  :: DataHandle
1797   CHARACTER*(*) :: Element
1798   CHARACTER*(*) :: VarName 
1799   real ,          INTENT(IN)  :: Data(*)
1800   INTEGER ,       INTENT(IN)  :: Count
1801   INTEGER ,       INTENT(OUT) :: Status
1802   CHARACTER(len=1000) :: tmpstr(1000)
1803   INTEGER             :: idx
1805   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_ti_real')
1807   if (committed(DataHandle)) then
1809      do idx = 1,Count
1810         write(tmpstr(idx),'(G17.10)')Data(idx)
1811      enddo
1813      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, Count, Status)
1815   endif
1817   RETURN
1818 END SUBROUTINE ext_gr1_put_var_ti_real 
1820 !*****************************************************************************
1822 SUBROUTINE ext_gr1_put_var_ti_double ( DataHandle,Element,  Varname, Data, &
1823      Count,  Status )
1824   USE gr1_data_info
1825   IMPLICIT NONE
1826 #include "wrf_status_codes.h"
1827   INTEGER ,       INTENT(IN)  :: DataHandle
1828   CHARACTER*(*) , INTENT(IN)  :: Element
1829   CHARACTER*(*) , INTENT(IN)  :: VarName
1830   real*8 ,            INTENT(IN) :: Data(*)
1831   INTEGER ,       INTENT(IN)  :: Count
1832   INTEGER ,       INTENT(OUT) :: Status
1833   CHARACTER(len=1000) :: tmpstr(1000)
1834   INTEGER             :: idx
1836   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_ti_double')
1838   if (committed(DataHandle)) then
1840      do idx = 1,Count
1841         write(tmpstr(idx),'(G17.10)')Data(idx)
1842      enddo
1843      
1844      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, Count, Status)
1845   endif
1847   RETURN
1848 END SUBROUTINE ext_gr1_put_var_ti_double
1850 !*****************************************************************************
1852 SUBROUTINE ext_gr1_put_var_ti_real8 ( DataHandle,Element,  Varname, Data, &
1853      Count,  Status )
1855   USE gr1_data_info
1856   IMPLICIT NONE
1857 #include "wrf_status_codes.h"
1858   INTEGER ,       INTENT(IN)  :: DataHandle
1859   CHARACTER*(*) :: Element
1860   CHARACTER*(*) :: VarName 
1861   real*8 ,        INTENT(IN)  :: Data(*)
1862   INTEGER ,       INTENT(IN)  :: Count
1863   INTEGER ,       INTENT(OUT) :: Status
1864   CHARACTER(len=1000) :: tmpstr(1000)
1865   INTEGER             :: idx
1867   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_ti_real8')
1869   if (committed(DataHandle)) then
1871      do idx = 1,Count
1872         write(tmpstr(idx),'(G17.10)')Data(idx)
1873      enddo
1874      
1875      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, Count, Status)
1876   endif
1878   RETURN
1879 END SUBROUTINE ext_gr1_put_var_ti_real8 
1881 !*****************************************************************************
1883 SUBROUTINE ext_gr1_put_var_ti_integer ( DataHandle,Element,  Varname, Data, &
1884      Count,  Status )
1886   USE gr1_data_info
1887   IMPLICIT NONE
1888 #include "wrf_status_codes.h"
1889   INTEGER ,       INTENT(IN)  :: DataHandle
1890   CHARACTER*(*) :: Element
1891   CHARACTER*(*) :: VarName 
1892   integer ,       INTENT(IN)  :: Data(*)
1893   INTEGER ,       INTENT(IN)  :: Count
1894   INTEGER ,       INTENT(OUT) :: Status
1895   CHARACTER(len=1000) :: tmpstr(1000)
1896   INTEGER             :: idx
1898   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_ti_integer')
1900   if (committed(DataHandle)) then
1902      do idx = 1,Count
1903         write(tmpstr(idx),'(G17.10)')Data(idx)
1904      enddo
1905      
1906      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, Count, Status)
1907   endif
1909   RETURN
1910 END SUBROUTINE ext_gr1_put_var_ti_integer 
1912 !*****************************************************************************
1914 SUBROUTINE ext_gr1_put_var_ti_logical ( DataHandle,Element,  Varname, Data, &
1915      Count,  Status )
1917   USE gr1_data_info
1918   IMPLICIT NONE
1919 #include "wrf_status_codes.h"
1920   INTEGER ,       INTENT(IN)  :: DataHandle
1921   CHARACTER*(*) :: Element
1922   CHARACTER*(*) :: VarName 
1923   logical ,       INTENT(IN)  :: Data(*)
1924   INTEGER ,       INTENT(IN)  :: Count
1925   INTEGER ,       INTENT(OUT) :: Status
1926   CHARACTER(len=1000) :: tmpstr(1000)
1927   INTEGER             :: idx
1929   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_ti_logical')
1931   if (committed(DataHandle)) then
1933      do idx = 1,Count
1934         write(tmpstr(idx),'(G17.10)')Data(idx)
1935      enddo
1936      
1937      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, Count, Status)
1939   endif
1941 RETURN
1942 END SUBROUTINE ext_gr1_put_var_ti_logical 
1944 !*****************************************************************************
1946 SUBROUTINE ext_gr1_put_var_ti_char ( DataHandle,Element,  Varname, Data,  &
1947      Status )
1949   USE gr1_data_info
1950   IMPLICIT NONE
1951 #include "wrf_status_codes.h"
1952   INTEGER ,       INTENT(IN)  :: DataHandle
1953   CHARACTER(len=*) :: Element
1954   CHARACTER(len=*) :: VarName 
1955   CHARACTER(len=*) :: Data
1956   INTEGER ,       INTENT(OUT) :: Status
1957   REAL dummy
1958   INTEGER                     :: Count
1959   CHARACTER(len=1000) :: tmpstr(1)
1960   INTEGER             :: idx
1962   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_ti_char')
1964   if (committed(DataHandle)) then
1966      write(tmpstr(1),*)trim(Data)
1968      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, 1, Status)
1970   endif
1972   RETURN
1973 END SUBROUTINE ext_gr1_put_var_ti_char 
1975 !******************************************************************************
1976 !* End of put_var_ti_* routines
1977 !******************************************************************************
1979 !******************************************************************************
1980 !* Start of get_var_td_* routines
1981 !******************************************************************************
1983 SUBROUTINE ext_gr1_get_var_td_double ( DataHandle,Element,  DateStr, &
1984      Varname, Data, Count, Outcount, Status )
1985   USE gr1_data_info
1986   IMPLICIT NONE
1987 #include "wrf_status_codes.h"
1988   INTEGER ,       INTENT(IN)  :: DataHandle
1989   CHARACTER*(*) , INTENT(IN)  :: Element
1990   CHARACTER*(*) , INTENT(IN)  :: DateStr
1991   CHARACTER*(*) , INTENT(IN)  :: VarName
1992   real*8 ,            INTENT(OUT) :: Data(*)
1993   INTEGER ,       INTENT(IN)  :: Count
1994   INTEGER ,       INTENT(OUT)  :: OutCount
1995   INTEGER ,       INTENT(OUT) :: Status
1996   INTEGER          :: idx
1997   INTEGER          :: stat
1998   CHARACTER*(1000) :: VALUE
2000   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_td_double')
2002   Status = WRF_NO_ERR
2003   
2004   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:),TRIM(Element),DateStr,&
2005        Varname,Value,stat)
2006   if (stat /= 0) then
2007      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2008      Status = WRF_WARN_VAR_NF
2009      RETURN
2010   endif
2012   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2013   if (stat .ne. 0) then
2014      CALL wrf_message("Reading data from"//Value//"failed")
2015      Status = WRF_WARN_COUNT_TOO_LONG
2016      RETURN
2017   endif
2018   Outcount = idx
2020 RETURN
2021 END SUBROUTINE ext_gr1_get_var_td_double
2023 !*****************************************************************************
2025 SUBROUTINE ext_gr1_get_var_td_real ( DataHandle,Element,  DateStr,Varname, &
2026      Data, Count, Outcount, Status )
2028   USE gr1_data_info
2029   IMPLICIT NONE
2030 #include "wrf_status_codes.h"
2031   INTEGER ,       INTENT(IN)  :: DataHandle
2032   CHARACTER*(*) :: Element
2033   CHARACTER*(*) :: DateStr
2034   CHARACTER*(*) :: VarName 
2035   real ,          INTENT(OUT) :: Data(*)
2036   INTEGER ,       INTENT(IN)  :: Count
2037   INTEGER ,       INTENT(OUT) :: OutCount
2038   INTEGER ,       INTENT(OUT) :: Status
2039   INTEGER          :: idx
2040   INTEGER          :: stat
2041   CHARACTER*(1000) :: VALUE
2043   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_td_real')
2045   Status = WRF_NO_ERR
2046   
2047   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), DateStr, &
2048        Varname, Value, stat)
2049   if (stat /= 0) then
2050      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2051      Status = WRF_WARN_VAR_NF
2052      RETURN
2053   endif
2055   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2056   if (stat .ne. 0) then
2057      CALL wrf_message("Reading data from"//Value//"failed")
2058      Status = WRF_WARN_COUNT_TOO_LONG
2059      RETURN
2060   endif
2061   Outcount = idx
2063   RETURN
2064 END SUBROUTINE ext_gr1_get_var_td_real 
2066 !*****************************************************************************
2068 SUBROUTINE ext_gr1_get_var_td_real8 ( DataHandle,Element,  DateStr,Varname, &
2069      Data, Count, Outcount, Status )
2071   USE gr1_data_info
2072   IMPLICIT NONE
2073 #include "wrf_status_codes.h"
2074   INTEGER ,       INTENT(IN)  :: DataHandle
2075   CHARACTER*(*) :: Element
2076   CHARACTER*(*) :: DateStr
2077   CHARACTER*(*) :: VarName 
2078   real*8 ,        INTENT(OUT) :: Data(*)
2079   INTEGER ,       INTENT(IN)  :: Count
2080   INTEGER ,       INTENT(OUT) :: OutCount
2081   INTEGER ,       INTENT(OUT) :: Status
2082   INTEGER          :: idx
2083   INTEGER          :: stat
2084   CHARACTER*(1000) :: VALUE
2086   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_td_real8')
2088   Status = WRF_NO_ERR
2089   
2090   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:),TRIM(Element),DateStr,&
2091        Varname,Value,stat)
2092   if (stat /= 0) then
2093      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2094      Status = WRF_WARN_VAR_NF
2095      RETURN
2096   endif
2098   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2099   if (stat .ne. 0) then
2100      CALL wrf_message("Reading data from"//Value//"failed")
2101      Status = WRF_WARN_COUNT_TOO_LONG
2102      RETURN
2103   endif
2104   Outcount = idx
2106   RETURN
2107 END SUBROUTINE ext_gr1_get_var_td_real8 
2109 !*****************************************************************************
2111 SUBROUTINE ext_gr1_get_var_td_integer ( DataHandle,Element,  DateStr,Varname, &
2112      Data, Count, Outcount, Status )
2114   USE gr1_data_info
2115   IMPLICIT NONE
2116 #include "wrf_status_codes.h"
2117   INTEGER ,       INTENT(IN)  :: DataHandle
2118   CHARACTER*(*) :: Element
2119   CHARACTER*(*) :: DateStr
2120   CHARACTER*(*) :: VarName 
2121   integer ,       INTENT(OUT) :: Data(*)
2122   INTEGER ,       INTENT(IN)  :: Count
2123   INTEGER ,       INTENT(OUT) :: OutCount
2124   INTEGER ,       INTENT(OUT) :: Status
2125   INTEGER          :: idx
2126   INTEGER          :: stat
2127   CHARACTER*(1000) :: VALUE
2129   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_td_integer')
2131   Status = WRF_NO_ERR
2132   
2133   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), DateStr, &
2134        Varname, Value,stat)
2135   if (stat /= 0) then
2136      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2137      Status = WRF_WARN_VAR_NF
2138      RETURN
2139   endif
2141   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2142   if (stat .ne. 0) then
2143      CALL wrf_message("Reading data from"//Value//"failed")
2144      Status = WRF_WARN_COUNT_TOO_LONG
2145      RETURN
2146   endif
2147   Outcount = idx
2149   RETURN
2150 END SUBROUTINE ext_gr1_get_var_td_integer 
2152 !*****************************************************************************
2154 SUBROUTINE ext_gr1_get_var_td_logical ( DataHandle,Element,  DateStr,Varname, &
2155      Data, Count, Outcount, Status )
2156   
2157   USE gr1_data_info
2158   IMPLICIT NONE
2159 #include "wrf_status_codes.h"
2160   INTEGER ,       INTENT(IN)  :: DataHandle
2161   CHARACTER*(*) :: Element
2162   CHARACTER*(*) :: DateStr
2163   CHARACTER*(*) :: VarName 
2164   logical ,       INTENT(OUT) :: Data(*)
2165   INTEGER ,       INTENT(IN)  :: Count
2166   INTEGER ,       INTENT(OUT) :: OutCount
2167   INTEGER ,       INTENT(OUT) :: Status
2168   INTEGER          :: idx
2169   INTEGER          :: stat
2170   CHARACTER*(1000) :: VALUE
2172   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_td_logical')
2174   Status = WRF_NO_ERR
2175   
2176   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), DateStr, &
2177        Varname, Value,stat)
2178   if (stat /= 0) then
2179      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2180      Status = WRF_WARN_VAR_NF
2181      RETURN
2182   endif
2184   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2185   if (stat .ne. 0) then
2186      CALL wrf_message("Reading data from"//Value//"failed")
2187      Status = WRF_WARN_COUNT_TOO_LONG
2188      RETURN
2189   endif
2190   Outcount = idx
2192   RETURN
2193 END SUBROUTINE ext_gr1_get_var_td_logical 
2195 !*****************************************************************************
2197 SUBROUTINE ext_gr1_get_var_td_char ( DataHandle,Element,  DateStr,Varname, &
2198      Data,  Status )
2200   USE gr1_data_info
2201   IMPLICIT NONE
2202 #include "wrf_status_codes.h"
2203   INTEGER ,       INTENT(IN)  :: DataHandle
2204   CHARACTER*(*) :: Element
2205   CHARACTER*(*) :: DateStr
2206   CHARACTER*(*) :: VarName 
2207   CHARACTER*(*) :: Data
2208   INTEGER ,       INTENT(OUT) :: Status
2209   INTEGER       :: stat
2211   Status = WRF_NO_ERR
2212   
2213   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_var_td_char')
2215   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), DateStr, &
2216        Varname, Data,stat)
2217   if (stat /= 0) then
2218      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2219      Status = WRF_WARN_VAR_NF
2220      RETURN
2221   endif
2223   RETURN
2224 END SUBROUTINE ext_gr1_get_var_td_char 
2226 !******************************************************************************
2227 !* End of get_var_td_* routines
2228 !******************************************************************************
2230 !******************************************************************************
2231 !* Start of put_var_td_* routines
2232 !******************************************************************************
2234 SUBROUTINE ext_gr1_put_var_td_double ( DataHandle, Element, DateStr, Varname, &
2235      Data, Count,  Status )
2236   USE gr1_data_info
2237   IMPLICIT NONE
2238 #include "wrf_status_codes.h"
2239   INTEGER ,       INTENT(IN)  :: DataHandle
2240   CHARACTER*(*) , INTENT(IN)  :: Element
2241   CHARACTER*(*) , INTENT(IN)  :: DateStr
2242   CHARACTER*(*) , INTENT(IN)  :: VarName
2243   real*8 ,            INTENT(IN) :: Data(*)
2244   INTEGER ,       INTENT(IN)  :: Count
2245   INTEGER ,       INTENT(OUT) :: Status
2246   CHARACTER(len=1000) :: tmpstr(1000)
2247   INTEGER             :: idx
2249   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_td_double')
2252   if (committed(DataHandle)) then
2254      do idx = 1,Count
2255         write(tmpstr(idx),'(G17.10)')Data(idx)
2256      enddo
2258      CALL gr1_build_string (td_output(DataHandle), &
2259           Varname//';'//DateStr//';'//Element, tmpstr, Count, Status)
2261   endif
2263 RETURN
2264 END SUBROUTINE ext_gr1_put_var_td_double
2266 !*****************************************************************************
2268 SUBROUTINE ext_gr1_put_var_td_integer ( DataHandle,Element,  DateStr, &
2269      Varname, Data, Count,  Status )
2271   USE gr1_data_info
2272   IMPLICIT NONE
2273 #include "wrf_status_codes.h"
2274   INTEGER ,       INTENT(IN)  :: DataHandle
2275   CHARACTER*(*) :: Element
2276   CHARACTER*(*) :: DateStr
2277   CHARACTER*(*) :: VarName 
2278   integer ,       INTENT(IN)  :: Data(*)
2279   INTEGER ,       INTENT(IN)  :: Count
2280   INTEGER ,       INTENT(OUT) :: Status
2281   CHARACTER(len=1000) :: tmpstr(1000)
2282   INTEGER             :: idx
2284   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_td_integer')
2286   if (committed(DataHandle)) then
2288      do idx = 1,Count
2289         write(tmpstr(idx),'(G17.10)')Data(idx)
2290      enddo
2291      
2292      CALL gr1_build_string (td_output(DataHandle), &
2293           Varname//';'//DateStr//';'//Element, tmpstr, Count, Status)
2295   endif
2297 RETURN
2298 END SUBROUTINE ext_gr1_put_var_td_integer 
2300 !*****************************************************************************
2302 SUBROUTINE ext_gr1_put_var_td_real ( DataHandle,Element,  DateStr,Varname, &
2303      Data, Count,  Status )
2305   USE gr1_data_info
2306   IMPLICIT NONE
2307 #include "wrf_status_codes.h"
2308   INTEGER ,       INTENT(IN)  :: DataHandle
2309   CHARACTER*(*) :: Element
2310   CHARACTER*(*) :: DateStr
2311   CHARACTER*(*) :: VarName 
2312   real ,          INTENT(IN)  :: Data(*)
2313   INTEGER ,       INTENT(IN)  :: Count
2314   INTEGER ,       INTENT(OUT) :: Status
2315   CHARACTER(len=1000) :: tmpstr(1000)
2316   INTEGER             :: idx
2318   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_td_real')
2320   if (committed(DataHandle)) then
2322      do idx = 1,Count
2323         write(tmpstr(idx),'(G17.10)')Data(idx)
2324      enddo
2325      
2326      CALL gr1_build_string (td_output(DataHandle), &
2327           Varname//';'//DateStr//';'//Element, tmpstr, Count, Status)
2329   endif
2331   RETURN
2332 END SUBROUTINE ext_gr1_put_var_td_real 
2334 !*****************************************************************************
2336 SUBROUTINE ext_gr1_put_var_td_real8 ( DataHandle,Element,  DateStr,Varname, &
2337      Data, Count,  Status )
2339   USE gr1_data_info
2340   IMPLICIT NONE
2341 #include "wrf_status_codes.h"
2342   INTEGER ,       INTENT(IN)  :: DataHandle
2343   CHARACTER*(*) :: Element
2344   CHARACTER*(*) :: DateStr
2345   CHARACTER*(*) :: VarName 
2346   real*8 ,        INTENT(IN)  :: Data(*)
2347   INTEGER ,       INTENT(IN)  :: Count
2348   INTEGER ,       INTENT(OUT) :: Status
2349   CHARACTER(len=1000) :: tmpstr(1000)
2350   INTEGER             :: idx
2352   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_td_real8')
2354   if (committed(DataHandle)) then
2355      do idx = 1,Count
2356         write(tmpstr(idx),'(G17.10)')Data(idx)
2357      enddo
2358      
2359      CALL gr1_build_string (td_output(DataHandle), &
2360           Varname//';'//DateStr//';'//Element, tmpstr, Count, Status)
2361   endif
2363   RETURN
2364 END SUBROUTINE ext_gr1_put_var_td_real8 
2366 !*****************************************************************************
2368 SUBROUTINE ext_gr1_put_var_td_logical ( DataHandle,Element,  DateStr, &
2369      Varname, Data, Count,  Status )
2371   USE gr1_data_info
2372   IMPLICIT NONE
2373 #include "wrf_status_codes.h"
2374   INTEGER ,       INTENT(IN)  :: DataHandle
2375   CHARACTER*(*) :: Element
2376   CHARACTER*(*) :: DateStr
2377   CHARACTER*(*) :: VarName 
2378   logical ,       INTENT(IN)  :: Data(*)
2379   INTEGER ,       INTENT(IN)  :: Count
2380   INTEGER ,       INTENT(OUT) :: Status
2381   CHARACTER(len=1000) :: tmpstr(1000)
2382   INTEGER             :: idx
2384   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_td_logical')
2386   if (committed(DataHandle)) then
2388      do idx = 1,Count
2389         write(tmpstr(idx),'(G17.10)')Data(idx)
2390      enddo
2392      CALL gr1_build_string (td_output(DataHandle), &
2393           Varname//';'//DateStr//';'//Element, tmpstr, Count, Status)
2395   endif
2397   RETURN
2398 END SUBROUTINE ext_gr1_put_var_td_logical 
2400 !*****************************************************************************
2402 SUBROUTINE ext_gr1_put_var_td_char ( DataHandle,Element,  DateStr,Varname, &
2403      Data,  Status )
2405   USE gr1_data_info
2406   IMPLICIT NONE
2407 #include "wrf_status_codes.h"
2408   INTEGER ,       INTENT(IN)  :: DataHandle
2409   CHARACTER*(*) :: Element
2410   CHARACTER*(*) :: DateStr
2411   CHARACTER*(*) :: VarName 
2412   CHARACTER*(*) :: Data
2413   INTEGER ,       INTENT(OUT) :: Status
2414   CHARACTER(len=1000) :: tmpstr(1)
2415   INTEGER             :: idx
2417   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_var_td_char')
2419   if (committed(DataHandle)) then
2421      write(tmpstr(idx),*)Data
2423      CALL gr1_build_string (td_output(DataHandle), &
2424           Varname//';'//DateStr//';'//Element, tmpstr, 1, Status)
2426   endif
2428   RETURN
2429 END SUBROUTINE ext_gr1_put_var_td_char 
2431 !******************************************************************************
2432 !* End of put_var_td_* routines
2433 !******************************************************************************
2436 !******************************************************************************
2437 !* Start of get_dom_ti_* routines
2438 !******************************************************************************
2440 SUBROUTINE ext_gr1_get_dom_ti_real ( DataHandle,Element,   Data, Count, &
2441      Outcount, Status )
2443   USE gr1_data_info
2444   IMPLICIT NONE
2445 #include "wrf_status_codes.h"
2446   INTEGER ,       INTENT(IN)  :: DataHandle
2447   CHARACTER*(*) :: Element
2448   real ,          INTENT(OUT) :: Data(*)
2449   INTEGER ,       INTENT(IN)  :: Count
2450   INTEGER ,       INTENT(OUT) :: Outcount
2451   INTEGER ,       INTENT(OUT) :: Status
2452   INTEGER          :: idx
2453   INTEGER          :: stat
2454   CHARACTER*(1000) :: VALUE
2456   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_ti_real')
2458   Status = WRF_NO_ERR
2459   
2460   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), "none", &
2461        "none", Value,stat)
2462   if (stat /= 0) then
2463      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2464      Status = WRF_WARN_VAR_NF
2465      RETURN
2466   endif
2468   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2469   if (stat .ne. 0) then
2470      CALL wrf_message("Reading data from"//Value//"failed")
2471      Status = WRF_WARN_COUNT_TOO_LONG
2472      RETURN
2473   endif
2474   Outcount = idx
2476   RETURN
2477 END SUBROUTINE ext_gr1_get_dom_ti_real 
2479 !*****************************************************************************
2481 SUBROUTINE ext_gr1_get_dom_ti_real8 ( DataHandle,Element,   Data, Count, &
2482      Outcount, Status )
2484   USE gr1_data_info
2485   IMPLICIT NONE
2486 #include "wrf_status_codes.h"
2487   INTEGER ,       INTENT(IN)  :: DataHandle
2488   CHARACTER*(*) :: Element
2489   real*8 ,        INTENT(OUT) :: Data(*)
2490   INTEGER ,       INTENT(IN)  :: Count
2491   INTEGER ,       INTENT(OUT) :: OutCount
2492   INTEGER ,       INTENT(OUT) :: Status
2493   INTEGER          :: idx
2494   INTEGER          :: stat
2495   CHARACTER*(1000) :: VALUE
2497   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_ti_real8')
2499   Status = WRF_NO_ERR
2500   
2501   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), "none", &
2502        "none", Value,stat)
2503   if (stat /= 0) then
2504      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2505      Status = WRF_WARN_VAR_NF
2506      RETURN
2507   endif
2509   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2510   if (stat .ne. 0) then
2511      CALL wrf_message("Reading data from"//Value//"failed")
2512      Status = WRF_WARN_COUNT_TOO_LONG
2513      RETURN
2514   endif
2515   Outcount = idx
2517   RETURN
2518 END SUBROUTINE ext_gr1_get_dom_ti_real8 
2520 !*****************************************************************************
2522 SUBROUTINE ext_gr1_get_dom_ti_integer ( DataHandle,Element,   Data, Count, &
2523      Outcount, Status )
2525   USE gr1_data_info
2526   IMPLICIT NONE
2527 #include "wrf_status_codes.h"
2528   INTEGER ,       INTENT(IN)  :: DataHandle
2529   CHARACTER*(*) :: Element
2530   integer ,       INTENT(OUT) :: Data(*)
2531   INTEGER ,       INTENT(IN)  :: Count
2532   INTEGER ,       INTENT(OUT) :: OutCount
2533   INTEGER ,       INTENT(OUT) :: Status
2534   INTEGER          :: idx
2535   INTEGER          :: stat
2536   CHARACTER*(1000) :: VALUE
2537   
2538   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_ti_integer Element: '//Element)
2540   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), "none", &
2541        "none", Value,stat)
2542   if (stat /= 0) then
2543      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2544      Status = WRF_WARN_VAR_NF
2545      RETURN
2546   endif
2548   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2549   if (stat .ne. 0) then
2550      CALL wrf_message("Reading data from"//Value//"failed")
2551      Status = WRF_WARN_COUNT_TOO_LONG
2552      RETURN
2553   endif
2554   Outcount = Count
2556   RETURN
2557 END SUBROUTINE ext_gr1_get_dom_ti_integer 
2559 !*****************************************************************************
2561 SUBROUTINE ext_gr1_get_dom_ti_logical ( DataHandle,Element,   Data, Count, &
2562      Outcount, Status )
2564   USE gr1_data_info
2565   IMPLICIT NONE
2566 #include "wrf_status_codes.h"
2567   INTEGER ,       INTENT(IN)  :: DataHandle
2568   CHARACTER*(*) :: Element
2569   logical ,       INTENT(OUT) :: Data(*)
2570   INTEGER ,       INTENT(IN)  :: Count
2571   INTEGER ,       INTENT(OUT) :: OutCount
2572   INTEGER ,       INTENT(OUT) :: Status
2573   INTEGER          :: idx
2574   INTEGER          :: stat
2575   CHARACTER*(1000) :: VALUE
2577   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_ti_logical')
2579   Status = WRF_NO_ERR
2580   
2581   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), "none", &
2582        "none", Value,stat)
2583   if (stat /= 0) then
2584      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2585      Status = WRF_WARN_VAR_NF
2586      RETURN
2587   endif
2589   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2590   if (stat .ne. 0) then
2591      CALL wrf_message("Reading data from"//Value//"failed")
2592      Status = WRF_WARN_COUNT_TOO_LONG
2593      RETURN
2594   endif
2595   Outcount = idx
2597   RETURN
2598 END SUBROUTINE ext_gr1_get_dom_ti_logical 
2600 !*****************************************************************************
2602 SUBROUTINE ext_gr1_get_dom_ti_char ( DataHandle,Element,   Data,  Status )
2604   USE gr1_data_info
2605   IMPLICIT NONE
2606 #include "wrf_status_codes.h"
2607   INTEGER ,       INTENT(IN)  :: DataHandle
2608   CHARACTER*(*) :: Element
2609   CHARACTER*(*) :: Data
2610   INTEGER ,       INTENT(OUT) :: Status
2611   INTEGER       :: stat
2612   INTEGER       :: endchar
2614   Status = WRF_NO_ERR
2615   
2616   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_ti_char')
2618   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), "none", &
2619        "none", Data, stat)
2620   if (stat /= 0) then
2621      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2622      Status = WRF_WARN_VAR_NF
2623      RETURN
2624   endif
2626   RETURN
2627 END SUBROUTINE ext_gr1_get_dom_ti_char 
2629 !*****************************************************************************
2631 SUBROUTINE ext_gr1_get_dom_ti_double ( DataHandle,Element,   Data, Count, &
2632      Outcount, Status )
2633   USE gr1_data_info
2634   IMPLICIT NONE
2635 #include "wrf_status_codes.h"
2636   INTEGER ,       INTENT(IN)  :: DataHandle
2637   CHARACTER*(*) , INTENT(IN)  :: Element
2638   real*8 ,            INTENT(OUT) :: Data(*)
2639   INTEGER ,       INTENT(IN)  :: Count
2640   INTEGER ,       INTENT(OUT)  :: OutCount
2641   INTEGER ,       INTENT(OUT) :: Status
2642   INTEGER          :: idx
2643   INTEGER          :: stat
2644   CHARACTER*(1000) :: VALUE
2646   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_ti_double')
2648   Status = WRF_NO_ERR
2649   
2650   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), "none", &
2651        "none", Value, stat)
2652   if (stat /= 0) then
2653      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2654      Status = WRF_WARN_VAR_NF
2655      RETURN
2656   endif
2658   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2659   if (stat .ne. 0) then
2660      CALL wrf_message("Reading data from"//Value//"failed")
2661      Status = WRF_WARN_COUNT_TOO_LONG
2662      RETURN
2663   endif
2664   Outcount = idx
2666 RETURN
2667 END SUBROUTINE ext_gr1_get_dom_ti_double
2669 !******************************************************************************
2670 !* End of get_dom_ti_* routines
2671 !******************************************************************************
2674 !******************************************************************************
2675 !* Start of put_dom_ti_* routines
2676 !******************************************************************************
2678 SUBROUTINE ext_gr1_put_dom_ti_real ( DataHandle,Element,   Data, Count,  &
2679      Status )
2681   USE gr1_data_info
2682   IMPLICIT NONE
2683 #include "wrf_status_codes.h"
2684   INTEGER ,       INTENT(IN)  :: DataHandle
2685   CHARACTER*(*) :: Element
2686   real ,          INTENT(IN)  :: Data(*)
2687   INTEGER ,       INTENT(IN)  :: Count
2688   INTEGER ,       INTENT(OUT) :: Status
2689   REAL dummy
2690   CHARACTER(len=1000) :: tmpstr(1000)
2691   character(len=2)    :: lf
2692   integer             :: idx
2694   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_ti_real')
2696   if (Element .eq. 'DX') then
2697      dx = Data(1)/1000.
2698   endif
2699   if (Element .eq. 'DY') then
2700      dy = Data(1)/1000.
2701   endif
2702   if (Element .eq. 'CEN_LAT') then
2703      center_lat = Data(1)
2704   endif
2705   if (Element .eq. 'CEN_LON') then
2706      center_lon = Data(1)
2707   endif  
2708   if (Element .eq. 'TRUELAT1') then
2709      truelat1 = Data(1)
2710   endif
2711   if (Element .eq. 'TRUELAT2') then
2712      truelat2 = Data(1)
2713   endif
2714   if (Element == 'STAND_LON') then
2715      proj_central_lon = Data(1)
2716   endif
2717   if (Element == 'DT') then
2718      timestep = Data(1)
2719   endif
2721   if (committed(DataHandle)) then
2723      do idx = 1,Count
2724         write(tmpstr(idx),'(G17.10)')Data(idx)
2725      enddo
2726      
2727      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, Count, Status)
2729   endif
2731   RETURN
2732 END SUBROUTINE ext_gr1_put_dom_ti_real 
2734 !*****************************************************************************
2736 SUBROUTINE ext_gr1_put_dom_ti_real8 ( DataHandle,Element,   Data, Count,  &
2737      Status )
2739   USE gr1_data_info
2740   IMPLICIT NONE
2741 #include "wrf_status_codes.h"
2742   INTEGER ,       INTENT(IN)  :: DataHandle
2743   CHARACTER*(*) :: Element
2744   real*8 ,        INTENT(IN)  :: Data(*)
2745   INTEGER ,       INTENT(IN)  :: Count
2746   INTEGER ,       INTENT(OUT) :: Status
2747   CHARACTER(len=1000) :: tmpstr(1000)
2748   INTEGER             :: idx
2750   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_ti_real8')
2752   if (committed(DataHandle)) then
2754      do idx = 1,Count
2755         write(tmpstr(idx),'(G17.10)')Data(idx)
2756      enddo
2757      
2758      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, Count, Status)
2760   endif
2762   RETURN
2763 END SUBROUTINE ext_gr1_put_dom_ti_real8 
2765 !*****************************************************************************
2767 SUBROUTINE ext_gr1_put_dom_ti_integer ( DataHandle,Element,   Data, Count,  &
2768      Status )
2770   USE gr1_data_info
2771   IMPLICIT NONE
2772 #include "wrf_status_codes.h"
2773   INTEGER ,       INTENT(IN)  :: DataHandle
2774   CHARACTER*(*) :: Element
2775   INTEGER ,       INTENT(IN)  :: Data(*)
2776   INTEGER ,       INTENT(IN)  :: Count
2777   INTEGER ,       INTENT(OUT) :: Status
2778   REAL dummy
2779   CHARACTER(len=1000) :: tmpstr(1000)
2780   INTEGER             :: idx
2783   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_ti_integer')
2785   if (Element == 'WEST-EAST_GRID_DIMENSION') then
2786      full_xsize = Data(1)
2787   else if (Element == 'SOUTH-NORTH_GRID_DIMENSION') then
2788      full_ysize = Data(1)
2789   else if (Element == 'MAP_PROJ') then
2790      projection = Data(1)
2791   else if (Element == 'WG_GRID_ID') then
2792      wg_grid_id = Data(1)
2793   endif
2795   if (committed(DataHandle)) then
2797      do idx = 1,Count
2798         write(tmpstr(idx),'(G17.10)')Data(idx)
2799      enddo
2800      
2801      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, Count, Status)
2803   endif
2805   call wrf_debug ( DEBUG , 'Leaving ext_gr1_put_dom_ti_integer')
2807   RETURN
2808 END SUBROUTINE ext_gr1_put_dom_ti_integer 
2810 !*****************************************************************************
2812 SUBROUTINE ext_gr1_put_dom_ti_logical ( DataHandle,Element,   Data, Count,  &
2813      Status )
2815   USE gr1_data_info
2816   IMPLICIT NONE
2817 #include "wrf_status_codes.h"
2818   INTEGER ,       INTENT(IN)  :: DataHandle
2819   CHARACTER*(*) :: Element
2820   logical ,       INTENT(IN)  :: Data(*)
2821   INTEGER ,       INTENT(IN)  :: Count
2822   INTEGER ,       INTENT(OUT) :: Status
2823   CHARACTER(len=1000) :: tmpstr(1000)
2824   INTEGER             :: idx
2826   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_ti_logical')
2828   if (committed(DataHandle)) then
2830      do idx = 1,Count
2831         write(tmpstr(idx),'(G17.10)')Data(idx)
2832      enddo
2833      
2834      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, Count, Status)
2836   endif
2838   RETURN
2839 END SUBROUTINE ext_gr1_put_dom_ti_logical 
2841 !*****************************************************************************
2843 SUBROUTINE ext_gr1_put_dom_ti_char ( DataHandle,Element,   Data,  &
2844      Status )
2846   USE gr1_data_info
2847   IMPLICIT NONE
2848 #include "wrf_status_codes.h"
2849   INTEGER ,       INTENT(IN)  :: DataHandle
2850   CHARACTER*(*) :: Element
2851   CHARACTER*(*),     INTENT(IN)  :: Data
2852   INTEGER ,       INTENT(OUT) :: Status
2853   REAL dummy
2854   CHARACTER(len=1000) :: tmpstr(1000)
2856   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_ti_char')
2858   if (Element .eq. 'START_DATE') then
2859      StartDate = Data
2860   endif
2862   if (committed(DataHandle)) then
2864      write(tmpstr(1),*)trim(Data)
2865      
2866      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, 1, Status)
2868   endif
2870   RETURN
2871 END SUBROUTINE ext_gr1_put_dom_ti_char
2873 !*****************************************************************************
2875 SUBROUTINE ext_gr1_put_dom_ti_double ( DataHandle,Element, Data, Count, &
2876      Status )
2877   USE gr1_data_info
2878   IMPLICIT NONE
2879 #include "wrf_status_codes.h"
2880   INTEGER ,       INTENT(IN)  :: DataHandle
2881   CHARACTER*(*) , INTENT(IN)  :: Element
2882   real*8 ,            INTENT(IN) :: Data(*)
2883   INTEGER ,       INTENT(IN)  :: Count
2884   INTEGER ,       INTENT(OUT) :: Status
2885   CHARACTER(len=1000) :: tmpstr(1000)
2886   INTEGER             :: idx
2888   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_ti_double')
2890   if (committed(DataHandle)) then
2892      do idx = 1,Count
2893         write(tmpstr(idx),'(G17.10)')Data(idx)
2894      enddo
2896      CALL gr1_build_string (ti_output(DataHandle), Element, tmpstr, Count, Status)
2898   endif
2899   
2900   RETURN
2901 END SUBROUTINE ext_gr1_put_dom_ti_double
2903 !******************************************************************************
2904 !* End of put_dom_ti_* routines
2905 !******************************************************************************
2908 !******************************************************************************
2909 !* Start of get_dom_td_* routines
2910 !******************************************************************************
2912 SUBROUTINE ext_gr1_get_dom_td_real ( DataHandle,Element, DateStr,  Data, &
2913      Count, Outcount, Status )
2915   USE gr1_data_info
2916   IMPLICIT NONE
2917 #include "wrf_status_codes.h"
2918   INTEGER ,       INTENT(IN)  :: DataHandle
2919   CHARACTER*(*) :: Element
2920   CHARACTER*(*) :: DateStr
2921   real ,          INTENT(OUT) :: Data(*)
2922   INTEGER ,       INTENT(IN)  :: Count
2923   INTEGER ,       INTENT(OUT) :: OutCount
2924   INTEGER ,       INTENT(OUT) :: Status
2925   INTEGER          :: idx
2926   INTEGER          :: stat
2927   CHARACTER*(1000) :: VALUE
2929   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_td_real')
2931   Status = WRF_NO_ERR
2932   
2933   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), DateStr, &
2934        "none", Value, stat)
2935   if (stat /= 0) then
2936      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2937      Status = WRF_WARN_VAR_NF
2938      RETURN
2939   endif
2941   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2942   if (stat .ne. 0) then
2943      CALL wrf_message("Reading data from"//Value//"failed")
2944      Status = WRF_WARN_COUNT_TOO_LONG
2945      RETURN
2946   endif
2947   Outcount = idx
2949   RETURN
2950 END SUBROUTINE ext_gr1_get_dom_td_real 
2952 !*****************************************************************************
2954 SUBROUTINE ext_gr1_get_dom_td_real8 ( DataHandle,Element, DateStr,  Data, &
2955      Count, Outcount, Status )
2957   USE gr1_data_info
2958   IMPLICIT NONE
2959 #include "wrf_status_codes.h"
2960   INTEGER ,       INTENT(IN)  :: DataHandle
2961   CHARACTER*(*) :: Element
2962   CHARACTER*(*) :: DateStr
2963   real*8 ,        INTENT(OUT) :: Data(*)
2964   INTEGER ,       INTENT(IN)  :: Count
2965   INTEGER ,       INTENT(OUT) :: OutCount
2966   INTEGER ,       INTENT(OUT) :: Status
2967   INTEGER          :: idx
2968   INTEGER          :: stat
2969   CHARACTER*(1000) :: VALUE
2971   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_td_real8')
2973   Status = WRF_NO_ERR
2974   
2975   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), DateStr, &
2976        "none", Value, stat)
2977   if (stat /= 0) then
2978      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
2979      Status = WRF_WARN_VAR_NF
2980      RETURN
2981   endif
2983   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
2984   if (stat .ne. 0) then
2985      CALL wrf_message("Reading data from"//Value//"failed")
2986      Status = WRF_WARN_COUNT_TOO_LONG
2987      RETURN
2988   endif
2989   Outcount = idx
2991   RETURN
2992 END SUBROUTINE ext_gr1_get_dom_td_real8 
2994 !*****************************************************************************
2996 SUBROUTINE ext_gr1_get_dom_td_integer ( DataHandle,Element, DateStr,  Data, &
2997      Count, Outcount, Status )
2999   USE gr1_data_info
3000   IMPLICIT NONE
3001 #include "wrf_status_codes.h"
3002   INTEGER ,       INTENT(IN)  :: DataHandle
3003   CHARACTER*(*) :: Element
3004   CHARACTER*(*) :: DateStr
3005   integer ,       INTENT(OUT) :: Data(*)
3006   INTEGER ,       INTENT(IN)  :: Count
3007   INTEGER ,       INTENT(OUT) :: OutCount
3008   INTEGER ,       INTENT(OUT) :: Status
3009   INTEGER          :: idx
3010   INTEGER          :: stat
3011   CHARACTER*(1000) :: VALUE
3013   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_td_integer')
3015   Status = WRF_NO_ERR
3016   
3017   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), DateStr, &
3018        "none", Value,stat)
3019   if (stat /= 0) then
3020      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
3021      Status = WRF_WARN_VAR_NF
3022      RETURN
3023   endif
3025   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
3026   if (stat .ne. 0) then
3027      CALL wrf_message("Reading data from"//Value//"failed")
3028      Status = WRF_WARN_COUNT_TOO_LONG
3029      RETURN
3030   endif
3031   Outcount = idx
3033   RETURN
3034 END SUBROUTINE ext_gr1_get_dom_td_integer 
3036 !*****************************************************************************
3038 SUBROUTINE ext_gr1_get_dom_td_logical ( DataHandle,Element, DateStr,  Data, &
3039      Count, Outcount, Status )
3041   USE gr1_data_info
3042   IMPLICIT NONE
3043 #include "wrf_status_codes.h"
3044   INTEGER ,       INTENT(IN)  :: DataHandle
3045   CHARACTER*(*) :: Element
3046   CHARACTER*(*) :: DateStr
3047   logical ,       INTENT(OUT) :: Data(*)
3048   INTEGER ,       INTENT(IN)  :: Count
3049   INTEGER ,       INTENT(OUT) :: OutCount
3050   INTEGER ,       INTENT(OUT) :: Status
3051   INTEGER          :: idx
3052   INTEGER          :: stat
3053   CHARACTER*(1000) :: VALUE
3055   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_td_logical')
3057   Status = WRF_NO_ERR
3058   
3059   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), DateStr, &
3060        "none", Value, stat)
3061   if (stat /= 0) then
3062      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
3063      Status = WRF_WARN_VAR_NF
3064      RETURN
3065   endif
3067   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
3068   if (stat .ne. 0) then
3069      CALL wrf_message("Reading data from"//Value//"failed")
3070      Status = WRF_WARN_COUNT_TOO_LONG
3071      RETURN
3072   endif
3073   Outcount = idx
3075   RETURN
3076 END SUBROUTINE ext_gr1_get_dom_td_logical 
3078 !*****************************************************************************
3080 SUBROUTINE ext_gr1_get_dom_td_char ( DataHandle,Element, DateStr,  Data,  &
3081      Status )
3083   USE gr1_data_info
3084   IMPLICIT NONE
3085 #include "wrf_status_codes.h"
3086   INTEGER ,       INTENT(IN)  :: DataHandle
3087   CHARACTER*(*) :: Element
3088   CHARACTER*(*) :: DateStr
3089   CHARACTER*(*) :: Data
3090   INTEGER ,       INTENT(OUT) :: Status
3091   INTEGER       :: stat
3093   Status = WRF_NO_ERR
3094   
3095   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_td_char')
3097   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), DateStr, &
3098        "none", Data, stat)
3099   if (stat /= 0) then
3100      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
3101      Status = WRF_WARN_VAR_NF
3102      RETURN
3103   endif
3105   RETURN
3106 END SUBROUTINE ext_gr1_get_dom_td_char 
3108 !*****************************************************************************
3110 SUBROUTINE ext_gr1_get_dom_td_double ( DataHandle,Element, DateStr,  Data, &
3111      Count, Outcount, Status )
3112   USE gr1_data_info
3113   IMPLICIT NONE
3114 #include "wrf_status_codes.h"
3115   INTEGER ,       INTENT(IN)  :: DataHandle
3116   CHARACTER*(*) , INTENT(IN)  :: Element
3117   CHARACTER*(*) , INTENT(IN)  :: DateStr
3118   real*8 ,            INTENT(OUT) :: Data(*)
3119   INTEGER ,       INTENT(IN)  :: Count
3120   INTEGER ,       INTENT(OUT)  :: OutCount
3121   INTEGER ,       INTENT(OUT) :: Status
3122   INTEGER          :: idx
3123   INTEGER          :: stat
3124   CHARACTER*(1000) :: VALUE
3126   call wrf_debug ( DEBUG , 'Entering ext_gr1_get_dom_td_double')
3128   Status = WRF_NO_ERR
3129   
3130   CALL GET_METADATA_VALUE(fileinfo(DataHandle)%fileindex(:), TRIM(Element), DateStr, &
3131        "none", Value, stat)
3132   if (stat /= 0) then
3133      CALL wrf_debug ( DEBUG , "GET_METADATA_VALUE failed for "//Element)
3134      Status = WRF_WARN_VAR_NF
3135      RETURN
3136   endif
3138   READ(Value,*,IOSTAT=stat)(Data(idx),idx=1,Count)
3139   if (stat .ne. 0) then
3140      CALL wrf_message("Reading data from"//Value//"failed")
3141      Status = WRF_WARN_COUNT_TOO_LONG
3142      RETURN
3143   endif
3144   Outcount = idx
3146 RETURN
3147 END SUBROUTINE ext_gr1_get_dom_td_double
3149 !******************************************************************************
3150 !* End of get_dom_td_* routines
3151 !******************************************************************************
3154 !******************************************************************************
3155 !* Start of put_dom_td_* routines
3156 !******************************************************************************
3159 SUBROUTINE ext_gr1_put_dom_td_real8 ( DataHandle,Element, DateStr,  Data, &
3160      Count,  Status )
3162   USE gr1_data_info
3163   IMPLICIT NONE
3164 #include "wrf_status_codes.h"
3165   INTEGER ,       INTENT(IN)  :: DataHandle
3166   CHARACTER*(*) :: Element
3167   CHARACTER*(*) :: DateStr
3168   real*8 ,        INTENT(IN)  :: Data(*)
3169   INTEGER ,       INTENT(IN)  :: Count
3170   INTEGER ,       INTENT(OUT) :: Status
3171   CHARACTER(len=1000) :: tmpstr(1000)
3172   INTEGER             :: idx
3174   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_td_real8')
3176   if (committed(DataHandle)) then
3178      do idx = 1,Count
3179         write(tmpstr(idx),'(G17.10)')Data(idx)
3180      enddo
3182      CALL gr1_build_string (td_output(DataHandle), DateStr//';'//Element, tmpstr, &
3183           Count, Status)
3185   endif
3187   RETURN
3188 END SUBROUTINE ext_gr1_put_dom_td_real8 
3190 !*****************************************************************************
3192 SUBROUTINE ext_gr1_put_dom_td_integer ( DataHandle,Element, DateStr,  Data, &
3193      Count,  Status )
3195   USE gr1_data_info
3196   IMPLICIT NONE
3197 #include "wrf_status_codes.h"
3198   INTEGER ,       INTENT(IN)  :: DataHandle
3199   CHARACTER*(*) :: Element
3200   CHARACTER*(*) :: DateStr
3201   integer ,       INTENT(IN)  :: Data(*)
3202   INTEGER ,       INTENT(IN)  :: Count
3203   INTEGER ,       INTENT(OUT) :: Status
3204   CHARACTER(len=1000) :: tmpstr(1000)
3205   INTEGER             :: idx
3207   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_td_integer')
3209   if (committed(DataHandle)) then
3211      do idx = 1,Count
3212         write(tmpstr(idx),'(G17.10)')Data(idx)
3213      enddo
3214      
3215      CALL gr1_build_string (td_output(DataHandle), DateStr//';'//Element, tmpstr, &
3216           Count, Status)
3218   endif
3220   RETURN
3221 END SUBROUTINE ext_gr1_put_dom_td_integer
3223 !*****************************************************************************
3225 SUBROUTINE ext_gr1_put_dom_td_logical ( DataHandle,Element, DateStr,  Data, &
3226      Count,  Status )
3228   USE gr1_data_info
3229   IMPLICIT NONE
3230 #include "wrf_status_codes.h"
3231   INTEGER ,       INTENT(IN)  :: DataHandle
3232   CHARACTER*(*) :: Element
3233   CHARACTER*(*) :: DateStr
3234   logical ,       INTENT(IN)  :: Data(*)
3235   INTEGER ,       INTENT(IN)  :: Count
3236   INTEGER ,       INTENT(OUT) :: Status
3237   CHARACTER(len=1000) :: tmpstr(1000)
3238   INTEGER             :: idx
3240   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_td_logical')
3242   if (committed(DataHandle)) then
3244      do idx = 1,Count
3245         write(tmpstr(idx),'(G17.10)')Data(idx)
3246      enddo
3247      
3248      CALL gr1_build_string (td_output(DataHandle), DateStr//';'//Element, tmpstr, &
3249           Count, Status)
3251   endif
3253   RETURN
3254 END SUBROUTINE ext_gr1_put_dom_td_logical
3256 !*****************************************************************************
3258 SUBROUTINE ext_gr1_put_dom_td_char ( DataHandle,Element, DateStr,  Data, &
3259      Status )
3261   USE gr1_data_info
3262   IMPLICIT NONE
3263 #include "wrf_status_codes.h"
3264   INTEGER ,       INTENT(IN)  :: DataHandle
3265   CHARACTER*(*) :: Element
3266   CHARACTER*(*) :: DateStr
3267   CHARACTER(len=*), INTENT(IN)  :: Data
3268   INTEGER ,       INTENT(OUT) :: Status
3269   CHARACTER(len=1000) :: tmpstr(1)
3271   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_td_char')
3273   if (committed(DataHandle)) then
3275      write(tmpstr(1),*)Data
3277      CALL gr1_build_string (td_output(DataHandle), DateStr//';'//Element, tmpstr, &
3278           1, Status)
3280   endif
3282   RETURN
3283 END SUBROUTINE ext_gr1_put_dom_td_char 
3285 !*****************************************************************************
3287 SUBROUTINE ext_gr1_put_dom_td_double ( DataHandle,Element, DateStr,  Data, &
3288      Count,  Status )
3289   USE gr1_data_info
3290   IMPLICIT NONE
3291 #include "wrf_status_codes.h"
3292   INTEGER ,       INTENT(IN)  :: DataHandle
3293   CHARACTER*(*) , INTENT(IN)  :: Element
3294   CHARACTER*(*) , INTENT(IN)  :: DateStr
3295   real*8 ,            INTENT(IN) :: Data(*)
3296   INTEGER ,       INTENT(IN)  :: Count
3297   INTEGER ,       INTENT(OUT) :: Status
3298   CHARACTER(len=1000) :: tmpstr(1000)
3299   INTEGER             :: idx
3301   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_td_double')
3303   if (committed(DataHandle)) then
3305      do idx = 1,Count
3306         write(tmpstr(idx),'(G17.10)')Data(idx)
3307      enddo
3309      CALL gr1_build_string (td_output(DataHandle), DateStr//';'//Element, tmpstr, &
3310           Count, Status)
3312   endif
3314 RETURN
3315 END SUBROUTINE ext_gr1_put_dom_td_double
3317 !*****************************************************************************
3319 SUBROUTINE ext_gr1_put_dom_td_real ( DataHandle,Element, DateStr,  Data, &
3320      Count,  Status )
3322   USE gr1_data_info
3323   IMPLICIT NONE
3324 #include "wrf_status_codes.h"
3325   INTEGER ,       INTENT(IN)  :: DataHandle
3326   CHARACTER*(*) :: Element
3327   CHARACTER*(*) :: DateStr
3328   real ,          INTENT(IN)  :: Data(*)
3329   INTEGER ,       INTENT(IN)  :: Count
3330   INTEGER ,       INTENT(OUT) :: Status
3331   CHARACTER(len=1000) :: tmpstr(1000)
3332   INTEGER             :: idx
3334   call wrf_debug ( DEBUG , 'Entering ext_gr1_put_dom_td_real')
3336   if (committed(DataHandle)) then
3338      do idx = 1,Count
3339         write(tmpstr(idx),'(G17.10)')Data(idx)
3340      enddo
3341      
3342      CALL gr1_build_string (td_output(DataHandle), DateStr//';'//Element, tmpstr, &
3343           Count, Status)
3345   endif
3347   RETURN
3348 END SUBROUTINE ext_gr1_put_dom_td_real 
3351 !******************************************************************************
3352 !* End of put_dom_td_* routines
3353 !******************************************************************************
3356 !*****************************************************************************
3358 SUBROUTINE gr1_build_string (string, Element, Value, Count, Status)
3360   IMPLICIT NONE
3361 #include "wrf_status_codes.h"
3363   CHARACTER (LEN=*) , INTENT(INOUT) :: string
3364   CHARACTER (LEN=*) , INTENT(IN)    :: Element
3365   CHARACTER (LEN=*) , INTENT(IN)    :: Value(*)
3366   INTEGER ,           INTENT(IN)    :: Count
3367   INTEGER ,           INTENT(OUT)   :: Status
3369   CHARACTER (LEN=2)                 :: lf
3370   INTEGER                           :: IDX
3372   lf=char(10)//' '
3373   if (len_trim(string) == 0) then
3374      string = lf//Element//' = '
3375   else
3376      string = trim(string)//lf//Element//' = '
3377   endif
3378   do idx = 1,Count
3379      if (idx > 1) then
3380         string = trim(string)//','
3381      endif
3382      string = trim(string)//' '//trim(adjustl(Value(idx)))
3383   enddo
3385   Status = WRF_NO_ERR
3387 END SUBROUTINE gr1_build_string
3389 !*****************************************************************************
3391 SUBROUTINE gr1_get_new_handle(DataHandle)
3392   USE gr1_data_info
3393   IMPLICIT NONE
3394   
3395   INTEGER ,       INTENT(OUT)  :: DataHandle
3396   INTEGER :: i
3398   DataHandle = -1
3399   do i=firstFileHandle, maxFileHandles
3400      if (.NOT. used(i)) then
3401         DataHandle = i
3402         used(i) = .true.
3403         exit
3404      endif
3405   enddo
3407   RETURN
3408 END SUBROUTINE gr1_get_new_handle
3411 !******************************************************************************
3414 SUBROUTINE gr1_get_levels(VarName, zidx, zsize, soil_layers, vert_stag, fraction, &
3415      vert_unit, level1, level2)
3417   use gr1_data_info
3418   IMPLICIT NONE
3420   integer :: zidx
3421   integer :: zsize
3422   logical :: soil_layers
3423   logical :: vert_stag
3424   logical :: fraction
3425   integer :: vert_unit
3426   integer :: level1
3427   integer :: level2
3428   character (LEN=*) :: VarName
3430   ! Setup vert_unit, and vertical levels in grib units
3432   if ((VarName .eq. 'LANDUSEF') .or. (VarName .eq. 'SOILCTOP') &
3433        .or. (VarName .eq. 'SOILCBOT')) then
3434      vert_unit = 109;
3435      level1 = zidx
3436      level2 = 0
3437   else if ((zsize .gt. 1) .and. (.not. soil_layers) .and. (.not. fraction)) &
3438        then
3439      vert_unit = 119;
3440      if (vert_stag) then
3441         level1 = (10000*full_eta(zidx)+0.5)
3442      else
3443         level1 = (10000*half_eta(zidx)+0.5)
3444      endif
3445      level2 = 0
3446   else
3447      ! Set the vertical coordinate and level for soil and 2D fields
3448      if (fraction) then
3449         vert_unit = 109
3450         level1 = zidx
3451         level2 = 0           
3452      else if (soil_layers) then
3453         vert_unit = 112
3454         level1 = 100*(soil_depth(zidx) - 0.5*soil_thickness(zidx))+0.5
3455         level2 = 100*(soil_depth(zidx) + 0.5*soil_thickness(zidx))+0.5
3456      else if (VarName .eq. 'mu') then
3457         vert_unit = 200
3458         level1 = 0
3459         level2 = 0
3460      else if ((VarName .eq. 'Q2') .or. (VarName .eq. 'TH2') .or. &
3461         (VarName .eq. 'T2')) then
3462         vert_unit = 105
3463         level1 = 2
3464         level2 = 0
3465      else if ((VarName .eq. 'Q10') .or. (VarName .eq. 'TH10') .or. &
3466           (VarName .eq. 'U10') .or. (VarName .eq. 'V10')) then
3467         vert_unit = 105
3468         level1 = 10
3469         level2 = 0
3470      else 
3471         vert_unit = 1
3472         level1 = 0
3473         level2 = 0
3474      endif
3475   endif
3477 end SUBROUTINE gr1_get_levels
3479 !*****************************************************************************
3482 SUBROUTINE gr1_fill_eta_levels(fileindex, FileFd, grib_tables, VarName, eta_levels)
3483   IMPLICIT NONE
3485   CHARACTER (len=*) :: fileindex
3486   INTEGER   :: FileFd
3487   CHARACTER (len=*) :: grib_tables
3488   character (len=*) :: VarName
3489   REAL,DIMENSION(*) :: eta_levels
3491   INTEGER   :: center, subcenter, parmtbl
3492   INTEGER   :: swapped
3493   INTEGER   :: leveltype
3494   INTEGER   :: idx
3495   INTEGER   :: parmid
3496   INTEGER   :: tablenum
3497   REAL      :: tmp
3498   INTEGER   :: numindices
3499   integer , DIMENSION(1000)   :: indices
3501   !
3502   ! Read the levels from the grib file
3503   !
3504   CALL GET_GRIB_PARAM(grib_tables, VarName, center, subcenter, parmtbl, &
3505        tablenum, parmid)
3507   if (parmid == -1) then
3508      call wrf_message ('Error getting grib parameter')
3509   endif
3511   leveltype = 119
3513   CALL GET_GRIB_INDICES(fileindex(:), center, subcenter, parmtbl, &
3514        parmid, "*", leveltype, &
3515        -HUGE(1), -HUGE(1), -HUGE(1), -HUGE(1), indices, numindices)
3517   do idx = 1,numindices
3518      CALL READ_GRIB(fileindex(:),FileFd,indices(idx),eta_levels(idx))
3519   enddo
3521   !
3522   ! Sort the levels--from highest (bottom) to lowest (top)
3523   !
3524   swapped = 1
3525   sortloop : do
3526      if (swapped /= 1) exit sortloop
3527      swapped = 0
3528      do idx=2, numindices
3529         !
3530         ! Remove duplicate levels, caused by multiple time periods in a 
3531         ! single file.
3532         !
3533         if (eta_levels(idx) == eta_levels(idx-1)) eta_levels(idx) = 0.0
3534         if (eta_levels(idx) > eta_levels(idx-1)) then
3535           tmp = eta_levels(idx)
3536           eta_levels(idx) = eta_levels(idx - 1)
3537           eta_levels(idx - 1) = tmp
3538           swapped = 1
3539         endif
3540      enddo
3541   enddo sortloop
3543 end subroutine gr1_fill_eta_levels