added README_changes.txt
[wrffire.git] / WPS / metgrid / src / write_met_module.F90
blobffdba350645e66058e9d98987f1cecc11f0a7b78
1 module write_met_module
3    use module_debug
4    use misc_definitions_module
5    use met_data_module
7    ! State variables?
8    integer :: output_unit
9    character (len=MAX_FILENAME_LEN) :: met_out_filename
11    contains
13    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
14    ! Name: write_met_init
15    !
16    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
17    subroutine write_met_init(fg_source, source_is_constant, datestr, istatus)
19       implicit none
20   
21       ! Arguments
22       integer, intent(out) :: istatus
23       logical, intent(in) :: source_is_constant
24       character (len=*), intent(in) :: fg_source
25       character (len=*), intent(in) :: datestr
26   
27       ! Local variables
28       integer :: io_status
29       logical :: is_used
31       istatus = 0
32     
33       !  1) BUILD FILENAME BASED ON TIME 
34       met_out_filename = ' '
35       if (.not. source_is_constant) then 
36          write(met_out_filename, '(a)') trim(fg_source)//':'//trim(datestr)
37       else
38          write(met_out_filename, '(a)') trim(fg_source)
39       end if
40   
41       !  2) OPEN FILE
42       do output_unit=10,100
43          inquire(unit=output_unit, opened=is_used)
44          if (.not. is_used) exit
45       end do 
46       call mprintf((output_unit > 100),ERROR,'In write_met_init(), couldn''t find an available Fortran unit.')
47       open(unit=output_unit, file=trim(met_out_filename), status='unknown', form='unformatted', iostat=io_status)
49       if (io_status > 0) istatus = 1
51       return
52   
54    end subroutine write_met_init
57    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
58    ! Name: write_next_met_field
59    !
60    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
61    subroutine write_next_met_field(fg_data, istatus)
63       implicit none
64   
65       ! Arguments
66       type (met_data), intent(in) :: fg_data
67       integer, intent(out) :: istatus
68   
69       ! Local variables
70       character (len=8) :: startloc
71       character (len=9) :: local_field
72   
73       istatus = 1
74   
75       !  1) WRITE FORMAT VERSION
76       write(unit=output_unit) fg_data % version
78       local_field = fg_data % field
79       if (local_field == 'GHT      ') local_field = 'HGT      '
81       ! PREGRID
82       if (fg_data % version == 3) then
84          ! Cylindrical equidistant
85          if (fg_data % iproj == PROJ_LATLON) then
86             write(unit=output_unit) fg_data % hdate, &
87                                     fg_data % xfcst, &
88                                     local_field,     &
89                                     fg_data % units, &
90                                     fg_data % desc,  &
91                                     fg_data % xlvl,  &
92                                     fg_data % nx,    &
93                                     fg_data % ny,    &
94                                     0
95             write(unit=output_unit) fg_data % startlat, &
96                                     fg_data % startlon, &
97                                     fg_data % deltalat, &
98                                     fg_data % deltalon
99      
100          ! Mercator
101          else if (fg_data % iproj == PROJ_MERC) then
102             write(unit=output_unit) fg_data % hdate, &
103                                     fg_data % xfcst, &
104                                     local_field,     &
105                                     fg_data % units, &
106                                     fg_data % desc,  &
107                                     fg_data % xlvl,  &
108                                     fg_data % nx,    &
109                                     fg_data % ny,    &
110                                     1
111             write(unit=output_unit) fg_data % startlat, &
112                                     fg_data % startlon, &
113                                     fg_data % dx,       &
114                                     fg_data % dy,       &
115                                     fg_data % truelat1
116      
117          ! Lambert conformal
118          else if (fg_data % iproj == PROJ_LC) then
119             write(unit=output_unit) fg_data % hdate, &
120                                     fg_data % xfcst, &
121                                     local_field,     &
122                                     fg_data % units, &
123                                     fg_data % desc,  &
124                                     fg_data % xlvl,  &
125                                     fg_data % nx,    &
126                                     fg_data % ny,    &
127                                     3
128             write(unit=output_unit) fg_data % startlat, &
129                                     fg_data % startlon, &
130                                     fg_data % dx,       &
131                                     fg_data % dy,       &
132                                     fg_data % xlonc,    &
133                                     fg_data % truelat1, &
134                                     fg_data % truelat2
135      
136          ! Polar stereographic
137          else if (fg_data % iproj == PROJ_PS) then
138             write(unit=output_unit) fg_data % hdate, &
139                                     fg_data % xfcst, &
140                                     local_field,     &
141                                     fg_data % units, &
142                                     fg_data % desc,  &
143                                     fg_data % xlvl,  &
144                                     fg_data % nx,    &
145                                     fg_data % ny,    &
146                                     5
147             write(unit=output_unit) fg_data % startlat, &
148                                     fg_data % startlon, &
149                                     fg_data % dx,       &
150                                     fg_data % dy,       &
151                                     fg_data % xlonc,    &
152                                     fg_data % truelat1
154          ! ?????????
155          else
156             call mprintf(.true.,ERROR,'Unrecognized projection code %i when reading from %s.', &
157                          i1=fg_data % iproj,s1=met_out_filename)
158      
159          end if
160      
161          write(unit=output_unit) fg_data % slab
162      
163          istatus = 0 
164     
165       ! GRIB_PREP
166       else if (fg_data % version == 4) then
168          if (fg_data % starti == 1.0 .and. fg_data % startj == 1.0) then
169             startloc='SWCORNER'
170          else
171             startloc='CENTER  '
172          end if
174          ! Cylindrical equidistant
175          if (fg_data % iproj == PROJ_LATLON) then
176             write(unit=output_unit) fg_data % hdate,      &
177                                     fg_data % xfcst,      &
178                                     fg_data % map_source, &
179                                     local_field,          &
180                                     fg_data % units,      &
181                                     fg_data % desc,       &
182                                     fg_data % xlvl,       &
183                                     fg_data % nx,         &
184                                     fg_data % ny,         &
185                                     0
186             write(unit=output_unit) startloc, &
187                                     fg_data % startlat, &
188                                     fg_data % startlon, &
189                                     fg_data % deltalat, &
190                                     fg_data % deltalon
192          ! Mercator
193          else if (fg_data % iproj == PROJ_MERC) then
194             write(unit=output_unit) fg_data % hdate,      &
195                                     fg_data % xfcst,      &
196                                     fg_data % map_source, &
197                                     local_field,          &
198                                     fg_data % units,      &
199                                     fg_data % desc,       &
200                                     fg_data % xlvl,       &
201                                     fg_data % nx,         &
202                                     fg_data % ny,         &
203                                     1
204             write(unit=output_unit) startloc, &
205                                     fg_data % startlat, &
206                                     fg_data % startlon, &
207                                     fg_data % dx,       &
208                                     fg_data % dy,       &
209                                     fg_data % truelat1
211          ! Lambert conformal
212          else if (fg_data % iproj == PROJ_LC) then
213             write(unit=output_unit) fg_data % hdate,      &
214                                     fg_data % xfcst,      &
215                                     fg_data % map_source, &
216                                     local_field,          &
217                                     fg_data % units,      &
218                                     fg_data % desc,       &
219                                     fg_data % xlvl,       &
220                                     fg_data % nx,         &
221                                     fg_data % ny,         &
222                                     3
223             write(unit=output_unit) startloc, &
224                                     fg_data % startlat, &
225                                     fg_data % startlon, &
226                                     fg_data % dx,       &
227                                     fg_data % dy,       &
228                                     fg_data % xlonc,    &
229                                     fg_data % truelat1, &
230                                     fg_data % truelat2
232          ! Polar stereographic
233          else if (fg_data % iproj == PROJ_PS) then
234             write(unit=output_unit) fg_data % hdate,      &
235                                     fg_data % xfcst,      &
236                                     fg_data % map_source, &
237                                     local_field,          &
238                                     fg_data % units,      &
239                                     fg_data % desc,       &
240                                     fg_data % xlvl,       &
241                                     fg_data % nx,         &
242                                     fg_data % ny,         &
243                                     5
244             write(unit=output_unit) startloc, &
245                                     fg_data % startlat, &
246                                     fg_data % startlon, &
247                                     fg_data % dx,       &
248                                     fg_data % dy,       &
249                                     fg_data % xlonc,    &
250                                     fg_data % truelat1
251      
252          ! ?????????
253          else
254             call mprintf(.true.,ERROR,'Unrecognized projection code %i when reading from %s.', &
255                          i1=fg_data % iproj,s1=met_out_filename)
256      
257          end if
258   
259          write(unit=output_unit) fg_data % slab
260       
261          istatus = 0
263       ! WPS
264       else if (fg_data % version == 5) then
266          if (fg_data % starti == 1.0 .and. fg_data % startj == 1.0) then
267             startloc='SWCORNER'
268          else
269             startloc='CENTER  '
270          end if
272          ! Cylindrical equidistant
273          if (fg_data % iproj == PROJ_LATLON) then
274             write(unit=output_unit) fg_data % hdate,      &
275                                     fg_data % xfcst,      &
276                                     fg_data % map_source, &
277                                     local_field,          &
278                                     fg_data % units,      &
279                                     fg_data % desc,       &
280                                     fg_data % xlvl,       &
281                                     fg_data % nx,         &
282                                     fg_data % ny,         &
283                                     0
284             write(unit=output_unit) startloc, &
285                                     fg_data % startlat, &
286                                     fg_data % startlon, &
287                                     fg_data % deltalat, &
288                                     fg_data % deltalon, &
289                                     fg_data % earth_radius
291          ! Mercator
292          else if (fg_data % iproj == PROJ_MERC) then
293             write(unit=output_unit) fg_data % hdate,      &
294                                     fg_data % xfcst,      &
295                                     fg_data % map_source, &
296                                     local_field,          &
297                                     fg_data % units,      &
298                                     fg_data % desc,       &
299                                     fg_data % xlvl,       &
300                                     fg_data % nx,         &
301                                     fg_data % ny,         &
302                                     1
303             write(unit=output_unit) startloc, &
304                                     fg_data % startlat, &
305                                     fg_data % startlon, &
306                                     fg_data % dx,       &
307                                     fg_data % dy,       &
308                                     fg_data % truelat1, &
309                                     fg_data % earth_radius
311          ! Lambert conformal
312          else if (fg_data % iproj == PROJ_LC) then
313             write(unit=output_unit) fg_data % hdate,      &
314                                     fg_data % xfcst,      &
315                                     fg_data % map_source, &
316                                     local_field,          &
317                                     fg_data % units,      &
318                                     fg_data % desc,       &
319                                     fg_data % xlvl,       &
320                                     fg_data % nx,         &
321                                     fg_data % ny,         &
322                                     3
323             write(unit=output_unit) startloc, &
324                                     fg_data % startlat, &
325                                     fg_data % startlon, &
326                                     fg_data % dx,       &
327                                     fg_data % dy,       &
328                                     fg_data % xlonc,    &
329                                     fg_data % truelat1, &
330                                     fg_data % truelat2, &
331                                     fg_data % earth_radius
333          ! Gaussian
334          else if (fg_data % iproj == PROJ_GAUSS) then
335             write(unit=output_unit) fg_data % hdate,      &
336                                     fg_data % xfcst,      &
337                                     fg_data % map_source, &
338                                     local_field,          &
339                                     fg_data % units,      &
340                                     fg_data % desc,       &
341                                     fg_data % xlvl,       &
342                                     fg_data % nx,         &
343                                     fg_data % ny,         &
344                                     4
345             write(unit=output_unit) startloc, &
346                                     fg_data % startlat, &
347                                     fg_data % startlon, &
348                                     fg_data % deltalat, &
349                                     fg_data % deltalon, &
350                                     fg_data % earth_radius
352          ! Polar stereographic
353          else if (fg_data % iproj == PROJ_PS) then
354             write(unit=output_unit) fg_data % hdate,      &
355                                     fg_data % xfcst,      &
356                                     fg_data % map_source, &
357                                     local_field,          &
358                                     fg_data % units,      &
359                                     fg_data % desc,       &
360                                     fg_data % xlvl,       &
361                                     fg_data % nx,         &
362                                     fg_data % ny,         &
363                                     5
364             write(unit=output_unit) startloc, &
365                                     fg_data % startlat, &
366                                     fg_data % startlon, &
367                                     fg_data % dx,       &
368                                     fg_data % dy,       &
369                                     fg_data % xlonc,    &
370                                     fg_data % truelat1, &
371                                     fg_data % earth_radius
372      
373          ! ?????????
374          else
375             call mprintf(.true.,ERROR,'Unrecognized projection code %i when reading from %s.', &
376                          i1=fg_data % iproj,s1=met_out_filename)
377      
378          end if
379   
380          write(unit=output_unit) fg_data % is_wind_grid_rel
382          write(unit=output_unit) fg_data % slab
383       
384          istatus = 0
386       else
387          call mprintf(.true.,ERROR,'Didn''t recognize format number %i.', i1=fg_data % version)
388       end if
389   
390       return
392    end subroutine write_next_met_field
395    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
396    ! Name: write_met_close
397    !
398    !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 
399    subroutine write_met_close()
401       implicit none
402   
403       close(unit=output_unit)
404       met_out_filename = 'UNINITIALIZED_FILENAME'
405   
406    end subroutine write_met_close
408 end module write_met_module