FV3 Bundle
station_data.F90
Go to the documentation of this file.
1 !***********************************************************************
2 !* GNU Lesser General Public License
3 !*
4 !* This file is part of the GFDL Flexible Modeling System (FMS).
5 !*
6 !* FMS is free software: you can redistribute it and/or modify it under
7 !* the terms of the GNU Lesser General Public License as published by
8 !* the Free Software Foundation, either version 3 of the License, or (at
9 !* your option) any later version.
10 !*
11 !* FMS is distributed in the hope that it will be useful, but WITHOUT
12 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
13 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
14 !* for more details.
15 !*
16 !* You should have received a copy of the GNU Lesser General Public
17 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
18 !***********************************************************************
19 
20 !> \author Giang Nong <Giang.Nong@gfdl.noaa.gov>
21 !!
22 !! \brief This module is used for outputing model results in a list
23 !! of stations (not gridded arrays). The user needs to supply
24 !! a list of stations with lat, lon values of each station.
25 !! Data at a single point (I,J) that is closest to each station will
26 !! be written to a file. No interpolation is made when a station
27 !! is between two or more grid points.<br>
28 !! In the output file, a 3D field will have a format of array(n1,n2) and
29 !! a 2D field is array(n1) where n1 is number of stations and n2 is number
30 !! of vertical levels or depths.
31 !!
32 !!
33 !! Here are some basic steps of how to use station_data_mod:<br>
34 !!
35 !! <ol>
36 !! <li> Call <tt>data_station_init</tt><br>
37 !! User needs to supply 2 tables: <tt>list_stations</tt> and <tt>station_data_table</tt> as follows:<br>
38 !! <table border="0">
39 !! <tr>
40 !! <td>Example of <tt>list_stations</tt>: &nbsp; &nbsp; &nbsp;</td>
41 !! <td>Example of <tt>station_data_table</tt>:</td>
42 !! </tr>
43 !! <tr>
44 !! <td>
45 !! <table border="0">
46 !! <tr>
47 !! <td>station_id &nbsp;</td>
48 !! <td>lat &nbsp;</td>
49 !! <td>lon</td>
50 !! </tr>
51 !! <tr>
52 !! <td>station_1</td>
53 !! <td>20.4</td>
54 !! <td>100.8</td>
55 !! </tr>
56 !! <tr>
57 !! <td>&nbsp;</td>
58 !! <td>&nbsp;</td>
59 !! <td>&nbsp;</td>
60 !! </tr>
61 !! <tr>
62 !! <td>&nbsp;</td>
63 !! <td>&nbsp;</td>
64 !! <td>&nbsp;</td>
65 !! </tr>
66 !! <tr>
67 !! <td>&nbsp;</td>
68 !! <td>&nbsp;</td>
69 !! <td>&nbsp;</td>
70 !! </tr>
71 !! <tr>
72 !! <td>&nbsp;</td>
73 !! <td>&nbsp;</td>
74 !! <td>&nbsp;</td>
75 !! </tr>
76 !! <tr>
77 !! <td>&nbsp;</td>
78 !! <td>&nbsp;</td>
79 !! <td>&nbsp;</td>
80 !! </tr>
81 !! <tr>
82 !! <td>&nbsp;</td>
83 !! <td>&nbsp;</td>
84 !! <td>&nbsp;</td>
85 !! </tr>
86 !! <tr>
87 !! <td>&nbsp;</td>
88 !! <td>&nbsp;</td>
89 !! <td>&nbsp;</td>
90 !! </tr>
91 !! <tr>
92 !! <td>&nbsp;</td>
93 !! <td>&nbsp;</td>
94 !! <td>&nbsp;</td>
95 !! </tr>
96 !! <tr>
97 !! <td>&nbsp;</td>
98 !! <td>&nbsp;</td>
99 !! <td>&nbsp;</td>
100 !! </tr>
101 !! </table>
102 !! </td>
103 !! <td>
104 !! <table border="0">
105 !! <tr>
106 !! <td>General descriptor:</td>
107 !! </tr>
108 !! <tr>
109 !! <td>&nbsp; &nbsp; Am2p14 station data</td>
110 !! </tr>
111 !! <tr>
112 !! <td>Start time (should be the same as model's initial time):</td>
113 !! </tr>
114 !! <tr>
115 !! <td>&nbsp; &nbsp; 19800101</td>
116 !! </tr>
117 !! <tr>
118 !! <td>File information:</td>
119 !! </tr>
120 !! <tr>
121 !! <td>
122 !! <table border="0">
123 !! <tr>
124 !! <td>filename &nbsp;</td>
125 !! <td>output_frequency &nbsp;</td>
126 !! <td>frequency_unit &nbsp;</td>
127 !! <td>time_axis_unit</td>
128 !! </tr>
129 !! <tr>
130 !! <td>"ocean_day"</td>
131 !! <td>1</td>
132 !! <td>"days"</td>
133 !! <td>"hours"</td>
134 !! </tr>
135 !! </table>
136 !! </td>
137 !! </tr>
138 !! <tr>
139 !! <td>Field information:</td>
140 !! </tr>
141 !! <tr>
142 !! <td>
143 !! <table border="0">
144 !! <tr>
145 !! <td>module &nbsp;</td>
146 !! <td>field_name &nbsp;</td>
147 !! <td>filename &nbsp;</td>
148 !! <td>time_method &nbsp;</td>
149 !! <td>pack</td>
150 !! </tr>
151 !! <tr>
152 !! <td>Ice_mod</td>
153 !! <td>temperature</td>
154 !! <td>ocean_day</td>
155 !! <td>.true.</td>
156 !! <td>2</td>
157 !! </tr>
158 !! <tr>
159 !! <td>Ice_mod</td>
160 !! <td>pressure</td>
161 !! <td>ocean_day</td>
162 !! <td>.false.</td>
163 !! <td>2</td>
164 !! </tr>
165 !! </table>
166 !! </td>
167 !! </tr>
168 !! </table>
169 !! </td>
170 !! </tr>
171 !! </table>
172 !! </li>
173 !! <li> Call <tt>register_station_field</tt> to register each field that needs to be written to a file, the call
174 !! <tt>register_station_field</tt> returns a field_id that will be used later in <tt>send_station_data</tt><br>
175 !! </li>
176 !! <li> Call <tt>send_station_data</tt> will send data at each station in the list
177 !! to a file <br>
178 !! </li>
179 !! <li> Finally, call <tt>station_data_end</tt> after the last time step.<br>
180 !! </li>
181 !! </ol>
182 !! Modules Included:
183 !! <table>
184 !! <tr>
185 !! <th>Module Name</th>
186 !! <th>Functions Included</th>
187 !! </tr>
188 !! <tr>
189 !! <td>axis_utils_mod</td>
190 !! <td>nearest_index</td>
191 !! </tr>
192 !! <tr>
193 !! <td>mpp_io_mod</td>
194 !! <td>mpp_open, MPP_RDONLY, MPP_ASCII, mpp_close, MPP_OVERRWR, MPP_NETCDF,
195 !! mpp_write_meta, MPP_SINGLE, mpp_write, fieldtype, mpp_flush</td>
196 !! </tr>
197 !! <tr>
198 !! <td>fms_mod</td>
199 !! <td>error_mesg, FATAL, WARNING, stdlog, write_version_number, mpp_pe, lowercase,
200 !! stdout, close_file, open_namelist_file, check_nml_error</td>
201 !! </tr>
202 !! <tr>
203 !! <td>mpp_mod</td>
204 !! <td>mpp_npes, mpp_sync, mpp_root_pe, mpp_send, mpp_recv, mpp_max,
205 !! mpp_get_current_pelist, input_nml_file, COMM_TAG_1, COMM_TAG_2,
206 !! COMM_TAG_3, COMM_TAG_4</td>
207 !! </tr>
208 !! <tr>
209 !! <td>mpp_domains_mod</td>
210 !! <td>domain2d, mpp_get_compute_domain</td>
211 !! </tr>
212 !! <tr>
213 !! <td>diag_axis_mod</td>
214 !! <td>diag_axis_init</td>
215 !! </tr>
216 !! <tr>
217 !! <td>diag_output_mod</td>
218 !! <td>write_axis_meta_data, write_field_meta_data, diag_fieldtype,
219 !! done_meta_data</td>
220 !! </tr>
221 !! <tr>
222 !! <td>diag_manager_mod</td>
223 !! <td>get_date_dif, DIAG_SECONDS, DIAG_MINUTES, DIAG_HOURS, DIAG_DAYS,
224 !! DIAG_MONTHS, DIAG_YEARS</td>
225 !! </tr>
226 !! <tr>
227 !! <td>diag_util_mod</td>
228 !! <td>diag_time_inc</td>
229 !! </tr>
230 !! <tr>
231 !! <td>time_manager_mod</td>
232 !! <td>operator(>), operator(>=), time_type, get_calendar_type, NO_CALENDAR,
233 !! set_time set_date, increment_date, increment_time</td>
234 !! </tr>
235 !! </table>
237 
238 
239 
240 use axis_utils_mod, only: nearest_index
241 use mpp_io_mod, only : mpp_open, mpp_rdonly, mpp_ascii, mpp_close,mpp_overwr,mpp_netcdf, &
242  mpp_write_meta, mpp_single, mpp_write, fieldtype,mpp_flush
243 use fms_mod, only : error_mesg, fatal, warning, stdlog, write_version_number,&
244  mpp_pe, lowercase, stdout, close_file, open_namelist_file, check_nml_error
245 use mpp_mod, only : mpp_npes, mpp_sync, mpp_root_pe, mpp_send, mpp_recv, mpp_max, &
246  mpp_get_current_pelist, input_nml_file, &
247  comm_tag_1, comm_tag_2, comm_tag_3, comm_tag_4
249 use diag_axis_mod, only : diag_axis_init
251 use diag_manager_mod,only : get_date_dif, diag_seconds, diag_minutes, diag_hours, &
252  diag_days, diag_months, diag_years
253 use diag_util_mod,only : diag_time_inc
254 use time_manager_mod, only: operator(>),operator(>=),time_type,get_calendar_type,no_calendar,set_time, &
256 
257 !--------------------------------------------------------------------
258 
259 implicit none
260 private
261 
262 !--------------------------------------------------------------------
263 !----------- Version number for this module -------------------------
264 
265 ! Include variable "version" to be written to log file.
266 #include<file_version.h>
267 
268 !--------------------------------------------------------------------
269 !------------------------- Namelist -----------------------------
270 
271 integer, parameter :: max_fields_per_file = 150
272 integer, parameter :: max_files = 31
273 integer :: num_files = 0
274 integer :: num_stations = 0
275 integer :: max_stations = 20
276 integer :: max_output_fields = 100
277 integer :: num_output_fields = 0
278 real :: empty = 0.0
279 real :: missing = 1.e20
280 logical :: module_is_initialized = .false.
281 logical :: need_write_axis = .true.
284 character (len=10) :: time_unit_list(6) = (/'seconds ', 'minutes ', &
285  'hours ', 'days ', 'months ', 'years '/)
286 integer, parameter :: every_time = 0
287 integer, parameter :: end_of_run = -1
288 
289 character(len=256) :: global_descriptor
290 character (len = 7) :: avg_name = 'average'
291 integer :: total_pe
292 integer :: lat_axis, lon_axis
293 integer, allocatable:: pelist(:)
294 character(len=32) :: pelist_name
296  character(len=128) :: name
297  integer :: output_freq
298  integer :: output_units
299  integer :: time_units
300  integer :: fields(max_fields_per_file)
301  integer :: num_fields
302  integer :: file_unit
303  integer :: time_axis_id, time_bounds_id
304  type(time_type) :: last_flush
305  type(fieldtype) :: f_avg_start, f_avg_end, f_avg_nitems, f_bounds
306 end type file_type
308  character(len=128) :: name
309  real :: lat, lon
310  integer :: id
311  integer :: global_i, global_j ! index of global grid
312  integer :: local_i, local_j ! index on the current PE
313  logical :: need_compute ! true if the station present in this PE
314 end type station_type
315 
317  integer :: output_file
318  integer :: num_station ! number of stations on this PE
319  integer, pointer :: station_id(:) =>null() ! id of station on this PE
320  character(len=128) :: output_name, module_name,long_name, units
321  logical :: time_average,time_max,time_min, time_ops, register
322  integer :: pack, axes(2), num_axes
323  character(len=8) :: time_method !could be: true, false, max, min, mean, ...
324  real, pointer :: buffer(:, :)=>null()
325  integer :: counter,nlevel
326  type(time_type) :: last_output, next_output
327  type(fieldtype) :: f_type
328 end type group_field_type
329 
331  real, pointer :: buffer(:,:)=>null()
332  integer :: counter
333 end type global_field_type
334 
335 
337 type(file_type),save :: files(max_files)
338 type(group_field_type),allocatable,save :: output_fields(:)
339 type(station_type),allocatable :: stations(:)
340 type(diag_fieldtype),save :: diag_field
342 
343 
344 !--------------------------------------------------------------------
345 !------------------------ Interfaces ----------------------------
346 
347 
348 !> \page register_station_field register_station_field Interface
349 !!
350 !! \brief register_station_field is similar to register_diag_field of diag_manager_mod.
351 !! All arguments are inputs that user needs to supply, some are optional. The
352 !! names of input args are self-describing. <br>
353 !! levels is absent for 2D fields. <br>
354 !! Note that pair (module_name, fieldname) must be unique in the
355 !! station_data_table of a fatal error will occur. <br>
356 !! A field id is returned from this call that will be used later in
357 !! send_station_data
359  module procedure register_station_field2d
360  module procedure register_station_field3d
361 end interface
362 
363 
364 !> \page send_station_data send_station_data Interface
365 !!
366 !! \brief Data should have the size of compute domain(isc:iec,jcs:jec)<br>
367 !! time is model's time<br>
368 !! field_id is returned from <tt>register_station_field</tt><br>
369 !! Only data at stations will be sent to root_pe which, in turn, sends
370 !! to output file
372  module procedure send_station_data_2d
373  module procedure send_station_data_3d
374 end interface
375 
376 !--------------------------------------------------------------------
377 
378 
379  contains
380 
381 
382 
383 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
384 !
385 ! Public Subroutines
386 !
387 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
388 
389 
390 
391 !> \brief Read in lat and lon of each station <br>
392 !! Create station_id based on lat, lon <br>
393 !! Read station_data_table, initialize output_fields and output files
394 subroutine station_data_init()
396 character(len=128) :: station_name
397 real :: lat, lon
398 integer :: iunit,nfiles,nfields,time_units,output_freq_units,j,station_id,io_status,logunit, ierr
399 logical :: init_verbose
400 character(len=128) :: record
401 type file_part_type
402  character(len=128) :: name
403  integer :: output_freq
404  character(len=10) :: output_freq_units
405  integer :: format ! must always be 1 for netcdf files
406  character(len=10) :: time_unit
407 end type file_part_type
408 type field_part_type
409  character(len=128) :: module_name,field_name,file_name
410  character(len=8) :: time_method
411  integer :: pack
412 end type field_part_type
413 
414 type(file_part_type) :: record_files
415 type(field_part_type) :: record_fields
416 
417 namelist /station_data_nml/ max_output_fields, max_stations,init_verbose
418 
419  if (module_is_initialized) return
420  init_verbose = .false.
421  total_pe = mpp_npes()
422  allocate(pelist(total_pe))
423  call mpp_get_current_pelist(pelist, pelist_name)
424 
425  !> read namelist
426 #ifdef INTERNAL_FILE_NML
427  read (input_nml_file, station_data_nml, iostat=io_status)
428  ierr = check_nml_error(io_status, 'station_data_nml')
429 #else
430  iunit = open_namelist_file()
431  ierr=1; do while (ierr /= 0)
432  read (iunit, nml=station_data_nml, iostat=io_status, end=10)
433  ierr = check_nml_error(io_status, 'station_data_nml')
434  enddo
435 10 call close_file (iunit)
436 
437 #endif
438  logunit = stdlog()
439  write(logunit, station_data_nml)
440 
442  !> read list of stations
443  if(init_verbose) then
444  logunit = stdout()
445  write(logunit, *) ' '
446  write(logunit, *) '****** Summary of STATION information from list_stations ********'
447  write(logunit, *) ' '
448  write(logunit, *) 'station name ', ' latitude ', ' longitude '
449  write(logunit, *) ' '
450  endif
451  call mpp_open(iunit, 'list_stations',form=mpp_ascii,action=mpp_rdonly)
452  do while (num_stations<max_stations)
453  read(iunit,'(a)',end=76,err=75) record
454  if (record(1:1) == '#') cycle
455  if(len_trim(record) < 1) cycle
456  read(record, *, end = 76, err = 75) station_name, lat, lon
457  station_id = get_station_id(lat, lon)
458  if(station_id > 0) then
459  stations(station_id)%name = station_name
460  else
461  call error_mesg('station_data_init','station DUPLICATED in file list_stations', fatal)
462  endif
463  logunit = stdout()
464  if( init_verbose.and. mpp_pe() == mpp_root_pe()) &
465  write(logunit,1)stations(station_id)%name,stations(station_id)%lat,stations(station_id)%lon
466 1 format(1x,a18, 1x,f8.2,4x,f8.2)
467 75 continue
468  enddo
469  call error_mesg('station_data_init','max_stations exceeded, increase it via namelist', fatal)
470 76 continue
471  call mpp_close (iunit)
472  logunit = stdout()
473  if(init_verbose) write(logunit, *)'*****************************************************************'
474 
475  !> read station_data table
476  call mpp_open(iunit, 'station_data_table',form=mpp_ascii,action=mpp_rdonly)
477  !> Read in the global file labeling string
478  read(iunit, *, end = 99, err=99) global_descriptor
479 
480  !> Read in the base date
481  read(iunit, *, end = 99, err = 99) base_year, base_month, base_day, &
483  if (get_calendar_type() /= no_calendar) then
486  else
487  !> If no calendar, ignore year and month
489  base_year = 0
490  base_month = 0
491  end if
492  nfiles=0
493  do while (nfiles <= max_files)
494  read(iunit,'(a)',end=86,err=85) record
495  if (record(1:1) == '#') cycle
496  read(record,*,err=85,end=85)record_files%name,record_files%output_freq, &
497  record_files%output_freq_units,record_files%format,record_files%time_unit
498  if(record_files%format /= 1) cycle !avoid reading field part
499  time_units = 0
500  output_freq_units = 0
501  do j = 1, size(time_unit_list(:))
502  if(record_files%time_unit == time_unit_list(j)) time_units = j
503  if(record_files%output_freq_units == time_unit_list(j)) output_freq_units = j
504  end do
505  if(time_units == 0) &
506  call error_mesg('station_data_init',' check time unit in station_data_table',fatal)
507  if(output_freq_units == 0) &
508  call error_mesg('station_data_init',', check output_freq in station_data_table',fatal)
509  call init_file(record_files%name,record_files%output_freq, output_freq_units,time_units)
510 85 continue
511  enddo
512  call error_mesg('station_data_init','max_files exceeded, increase max_files', fatal)
513 86 continue
514  rewind(iunit)
515  nfields=0
516  do while (nfields <= max_output_fields)
517  read(iunit,'(a)',end=94,err=93) record
518  if (record(1:1) == '#') cycle
519  read(record,*,end=93,err=93) record_fields
520  if (record_fields%pack .gt. 8 .or.record_fields%pack .lt. 1) cycle !avoid reading file part
521  nfields=nfields+1
522  call init_output_field(record_fields%module_name,record_fields%field_name, &
523  record_fields%file_name,record_fields%time_method,record_fields%pack)
524 93 continue
525  enddo
526  call error_mesg('station_data_init','max_output_fields exceeded, increase it via nml ', fatal)
527 94 continue
528  call close_file(iunit)
530  call write_version_number ("STATION_DATA_MOD", version)
531  module_is_initialized = .true.
532  return
533 99 continue
534  call error_mesg('station_data_init','error reading station_datatable',fatal)
535 end subroutine station_data_init
536 
537 
538 
539 
540 !> \brief check_duplicate_output_fields takes the data pairs
541 !! (output_name and output_file) and (module_name and output_name) and
542 !! ensures each pair is unique.
543 !!
544 !! \throw FATAL, "ERROR1 in station_data_table: module/field <module/name> duplicated"
545 !! \throw FATAL, "ERROR2 in station_data_table: module/field <module/name> duplicated"
547 ! pair(output_name and output_file) should be unique in data_station_table, ERROR1
548 ! pair(module_name and output_name) should be unique in data_station_table, ERROR2
549 integer :: i, j, tmp_file
550 character(len=128) :: tmp_name, tmp_module
551 
552 if(num_output_fields <= 1) return
553 do i = 1, num_output_fields-1
554  tmp_name = trim(output_fields(i)%output_name)
555  tmp_file = output_fields(i)%output_file
556  tmp_module = trim(output_fields(i)%module_name)
557  do j = i+1, num_output_fields
558  if((tmp_name == trim(output_fields(j)%output_name)).and. &
559  (tmp_file == output_fields(j)%output_file)) &
560  call error_mesg (' ERROR1 in station_data_table:', &
561  &' module/field '//tmp_module//'/'//tmp_name//' duplicated', fatal)
562  if((tmp_name == trim(output_fields(j)%output_name)).and. &
563  (tmp_module == trim(output_fields(j)%module_name))) &
564  call error_mesg (' ERROR2 in station_data_table:', &
565  &' module/field '//tmp_module//'/'//tmp_name//' duplicated', fatal)
566  enddo
567 enddo
568 end subroutine check_duplicate_output_fields
569 
570 
571 !> \brief get_station_id is passed the station's distinct lat and lon
572 !! to determine what the station's id is.
573 function get_station_id(lat,lon)
574  integer :: i
575  integer :: get_station_id !< Unique ID of station
576  real, intent(in):: lat !< Latitude of station
577  real, intent(in):: lon !< Longitude of station
578 ! each station should have distinct lat and lon
579  get_station_id = -1
580  do i = 1, num_stations
581  if(stations(i)%lat == lat .and. stations(i)%lon == lon) return
582  enddo
585  stations(num_stations)%lat = lat
586  stations(num_stations)%lon = lon
587  stations(num_stations)%need_compute = .false.
588  stations(num_stations)%global_i = -1; stations(num_stations)%global_j = -1
589  stations(num_stations)%local_i = -1 ; stations(num_stations)%local_j = -1
591 end function get_station_id
592 
593 
594 !> \brief init_file initializes files to write station data to. Data is formatted as
595 !! number of units of measurement at X output frequency over amount of time dt
596 !!
597 !! \throw FATAL, "init_file or max_files exceeded, increase max_files"
598 subroutine init_file(filename, output_freq, output_units, time_units)
599  character(len=*), intent(in) :: filename !< Filename to use when initializing file
600  integer, intent(in) :: output_freq !< Output frequency of substance
601  integer, intent(in) :: output_units !< Number of units of substance being recorded
602  integer, intent(in) :: time_units !< Length of time being recorded
603  character(len=128) :: time_units_str
604  real, dimension(1) :: tdata
605 
606  num_files = num_files + 1
607  if(num_files >= max_files) &
608  call error_mesg('station_data, init_file', ' max_files exceeded, incease max_files', fatal)
609  files(num_files)%name = trim(filename)
610  files(num_files)%output_freq = output_freq
611  files(num_files)%output_units = output_units
612  files(num_files)%time_units = time_units
613  files(num_files)%num_fields = 0
614  files(num_files)%last_flush = base_time
615  files(num_files)%file_unit = -1
616  !> register axis_id and time boundaries id
617  write(time_units_str, 11) trim(time_unit_list(files(num_files)%time_units)), base_year, &
619 11 format(a, ' since ', i4.4, '-', i2.2, '-', i2.2, ' ', i2.2, ':', i2.2, ':', i2.2)
620  files(num_files)%time_axis_id = diag_axis_init('Time', tdata, time_units_str, 'T', &
621  'Time' , set_name=trim(filename))
622  files(num_files)%time_bounds_id = diag_axis_init('nv',(/1.,2./),'none','N','vertex number',&
623  set_name=trim(filename))
624 end subroutine init_file
625 
626 
627 !> \brief init_output_field initializes output_field attributes
628 !!
629 !! \throw FATAL, "max_output_fields exceeded, increase it via nml"
630 !! \throw FATAL, "file <filename.ext> is NOT found in station_data_table"
631 !! \throw FATAL, "max_fields_per_file exceeded"
632 !! \throw FATAL, "time_method MAX is not supported"
633 !! \throw FATAL, "time_method MIN is not supported"
634 !! \throw FATAL, "error in time_method of field field_name"
635 subroutine init_output_field(module_name,field_name,file_name,time_method,pack)
636  character(len=*), intent(in) :: module_name, field_name, file_name
637  character(len=*), intent(in) :: time_method
638  integer, intent(in) :: pack
639  integer :: out_num, file_num,num_fields, method_selected, l1
640  character(len=8) :: t_method
641  !> Get a number for this output field
644  call error_mesg('station_data', 'max_output_fields exceeded, increase it via nml', fatal)
645  out_num = num_output_fields
646  file_num = find_file(file_name)
647  if(file_num < 0) &
648  call error_mesg('station_data,init_output_field', 'file '//trim(file_name) &
649  //' is NOT found in station_data_table', fatal)
650  !> Insert this field into list of fields of this file
651  files(file_num)%num_fields = files(file_num)%num_fields + 1
652  if(files(file_num)%num_fields > max_fields_per_file) &
653  call error_mesg('station_data, init_output_field', 'max_fields_per_file exceeded ', fatal)
654  num_fields = files(file_num)%num_fields
655  files(file_num)%fields(num_fields) = out_num
656  output_fields(out_num)%output_name = trim(field_name)
657  output_fields(out_num)%module_name = trim(module_name)
658  output_fields(out_num)%counter = 0
659  output_fields(out_num)%output_file = file_num
660  output_fields(out_num)%pack = pack
661  output_fields(out_num)%time_average = .false.
662  output_fields(out_num)%time_min = .false.
663  output_fields(out_num)%time_max = .false.
664  output_fields(out_num)%time_ops = .false.
665  output_fields(out_num)%register = .false.
666  t_method = lowercase(time_method)
667  select case (trim(t_method))
668  case('.true.')
669  output_fields(out_num)%time_average = .true.
670  output_fields(out_num)%time_method = 'mean'
671  case('mean')
672  output_fields(out_num)%time_average = .true.
673  output_fields(out_num)%time_method = 'mean'
674  case('average')
675  output_fields(out_num)%time_average = .true.
676  output_fields(out_num)%time_method = 'mean'
677  case('avg')
678  output_fields(out_num)%time_average = .true.
679  output_fields(out_num)%time_method = 'mean'
680  case('.false.')
681  output_fields(out_num)%time_average = .false.
682  output_fields(out_num)%time_method = 'point'
683  case ('max')
684  call error_mesg('station_data, init_output_field','time_method MAX is not supported',&
685  fatal)
686  output_fields(out_num)%time_max = .true.
687  output_fields(out_num)%time_method = 'max'
688  l1 = len_trim(output_fields(out_num)%output_name)
689  if(output_fields(out_num)%output_name(l1-2:l1) /= 'max') &
690  output_fields(out_num)%output_name = trim(field_name)//'_max'
691  case ('min')
692  call error_mesg('station_data, init_output_field','time_method MIN is not supported',&
693  fatal)
694  output_fields(out_num)%time_min = .true.
695  output_fields(out_num)%time_method = 'min'
696  l1 = len_trim(output_fields(out_num)%output_name)
697  if(output_fields(out_num)%output_name(l1-2:l1) /= 'min') &
698  output_fields(out_num)%output_name = trim(field_name)//'_min'
699  case default
700  call error_mesg('station_data, init_output_field', 'error in time_method of field '&
701  //trim(field_name), fatal)
702  end select
703  if (files(file_num)%output_freq == every_time) &
704  output_fields(out_num)%time_average = .false.
705  output_fields(out_num)%time_ops = output_fields(out_num)%time_min.or.output_fields(out_num)%time_max &
706  .or.output_fields(out_num)%time_average
707  output_fields(out_num)%time_method = trim(time_method)
708 end subroutine init_output_field
709 
710 !> \brief find_file finds the index of a requested file
711 !!
712 !! \param [in] <name> Name of the file
713 function find_file(name)
714 integer :: find_file
715 character(len=*), intent(in) :: name
716 integer :: i
717 
718 find_file = -1
719 do i = 1, num_files
720  if(trim(files(i)%name) == trim(name)) then
721  find_file = i
722  return
723  end if
724 end do
725 end function find_file
726 
727 
728 !> \brief register_station_field2d registers a new station field with user input.<br>
729 !! register_station_field2d is overloaded by register_station_field3d as a part
730 !! of the register_station_field interface.
731 function register_station_field2d (module_name,fieldname,glo_lat,glo_lon,init_time, &
732  domain,longname,units)
734  character(len=*), intent(in) :: module_name !< Module name of the station being registered
735  character(len=*), intent(in) :: fieldname !< Field name of the tation being registered
736  real,dimension(:), intent(in) :: glo_lat !< Global latitude location of the station
737  real,dimension(:), intent(in) :: glo_lon !< Global longitude location of the station
738  type(domain2d), intent(in) :: domain !< Optional 2d domain representation
739  type(time_type), intent(in) :: init_time !< Initial time of the station
740  character(len=*), optional, intent(in) :: longname !< Optional longname to identify station
741  character(len=*), optional, intent(in) :: units !< Optional units for station measurement data
742  real :: levels(1:1)
743 
744  levels = 0.
745  register_station_field2d = register_station_field3d(module_name,fieldname,glo_lat,glo_lon,&
746  levels,init_time,domain,longname,units)
747 end function register_station_field2d
748 
749 
750 !> \brief register_station_field3d registers a new station field with user input.<br>
751 !! register_station_field3d is overloaded by register_station_field2d as a part
752 !! of the register_station_field interface.
753 !!
754 !! \throw FATAL, "outside global latitude values"
755 !! \throw FATAL, "outside global longitude values"
756 !! \throw FATAL, "Error in global index of station"
757 !! \throw WARNING, "<module_name>/<field_name> NOT found in station_data table"
758 !! \throw FATAL, "Error in determining local_num_station"
759 function register_station_field3d (module_name,fieldname,glo_lat,glo_lon,levels,init_time, &
760  domain,longname,units)
762  !> write field meta data on ROOT PE only
763  !> allocate buffer
764  integer :: register_station_field3d
765  character(len=*), intent(in) :: module_name !< Module name of the station being registered
766  character(len=*), intent(in) :: fieldname !< Field name of the tation being registered
767  real,dimension(:), intent(in) :: glo_lat !< Global latitude location of the station (in X direction)
768  real,dimension(:), intent(in) :: glo_lon !< Global longitude location of the station (in Y direction)
769  real,dimension(:), intent(in) :: levels !< Global elevation location of the station (in Z direction)
770  type(domain2d), intent(in) :: domain !< Optional 2d domain representation
771  type(time_type), intent(in) :: init_time !< Initial time of the station
772  character(len=*), optional, intent(in) :: longname !< Optional longname to identify station
773  character(len=*), optional, intent(in) :: units !< Optional units for station measurement data
774  integer :: i,ii, nlat, nlon,nlevel, isc, iec, jsc, jec
775  character(len=128) :: error_msg
776  integer :: local_num_stations !< number of stations on this PE
777  integer :: out_num !< position of this field in array output_fields
778  integer :: file_num, freq, output_units, outunit
779  real, allocatable :: station_values(:), level_values(:)
780  character(len=128) :: longname2,units2
781 
782 
783  if(PRESENT(longname)) then
784  longname2 = longname
785  else
786  longname2 = fieldname
787  endif
788  if(PRESENT(units)) then
789  units2 = units
790  else
791  units2 = "none"
792  endif
793 
794  nlat = size(glo_lat); nlon = size(glo_lon); nlevel=size(levels)
795  allocate(station_values(num_stations), level_values(nlevel))
796  do i = 1, nlevel
797  level_values(i) = real(i)
798  enddo
799  !> determine global index of this field in all stations
800  outunit = stdout()
801  do i = 1,num_stations
802  station_values(i) = real(i)
803  if(stations(i)%lat<glo_lat(1) .or. stations(i)%lat>glo_lat(nlat)) then
804  write(error_msg,'(F9.3)') stations(i)%lat
805  write(outunit,*) 'Station with latitude '//trim(error_msg)//' outside global latitude values'
806  call error_mesg ('register_station_field', 'latitude out of range', fatal)
807  endif
808  if(stations(i)%lon<glo_lon(1) .or. stations(i)%lon>glo_lon(nlon)) then
809  write(error_msg,'(F9.3)') stations(i)%lon
810  write(outunit,*) 'Station with longitude '//trim(error_msg)//' outside global longitude values'
811  call error_mesg ('register_station_field', 'longitude out of range', fatal)
812  endif
813  stations(i)%global_i = nearest_index(stations(i)%lon, glo_lon)
814  stations(i)%global_j = nearest_index(stations(i)%lat, glo_lat)
815  if(stations(i)%global_i<0 .or. stations(i)%global_j<0) &
816  call error_mesg ('register_station_field', 'Error in global index of station',fatal)
817  enddo
818  !> determine local index of this field in all stations , local index starts from 1
819  call mpp_get_compute_domain(domain, isc,iec,jsc,jec)
820  local_num_stations = 0
821  do i = 1,num_stations
822  if(isc<=stations(i)%global_i .and. iec>= stations(i)%global_i .and. &
823  jsc<=stations(i)%global_j .and. jec>= stations(i)%global_j) then
824  stations(i)%need_compute = .true.
825  stations(i)%local_i = stations(i)%global_i - isc + 1
826  stations(i)%local_j = stations(i)%global_j - jsc + 1
827  local_num_stations = local_num_stations +1
828  endif
829  enddo
830  !> get the position of this field in the array output_fields
831  out_num = find_output_field(module_name, fieldname)
832  if(out_num < 0 .and. mpp_pe() == mpp_root_pe()) then
833  call error_mesg ('register_station_field', &
834  'module/field_name '//trim(module_name)//'/'//&
835  trim(fieldname)//' NOT found in station_data table', warning)
836  register_station_field3d = out_num
837  return
838  endif
839  if(local_num_stations>0) then
840  allocate(output_fields(out_num)%station_id(local_num_stations))
841  allocate(output_fields(out_num)%buffer(local_num_stations,nlevel))
842  output_fields(out_num)%buffer = empty
843  !> fill out list of available stations in this PE
844  ii=0
845  do i = 1,num_stations
846  if(stations(i)%need_compute) then
847  ii = ii+ 1
848  if(ii>local_num_stations) call error_mesg ('register_station_field', &
849  'error in determining local_num_station', fatal)
850  output_fields(out_num)%station_id(ii)=stations(i)%id
851  endif
852  enddo
853  endif
854  output_fields(out_num)%num_station = local_num_stations
855  if( mpp_pe() == mpp_root_pe()) then
856  allocate(global_field%buffer(num_stations,nlevel))
857  global_field%buffer = missing
858  endif
859  output_fields(out_num)%register = .true.
860  output_fields(out_num)%output_name = fieldname
861  file_num = output_fields(out_num)%output_file
862  output_fields(out_num)%last_output = init_time
863  freq = files(file_num)%output_freq
864  output_units = files(file_num)%output_units
865  output_fields(out_num)%next_output = diag_time_inc(init_time, freq, output_units)
866  register_station_field3d = out_num
867  output_fields(out_num)%long_name = longname2
868  output_fields(out_num)%units = units2
869  output_fields(out_num)%nlevel = nlevel
870  !> deal with axes
871 
872  output_fields(out_num)%axes(1) = diag_axis_init('Stations',station_values,'station number', 'X')
873  if(nlevel == 1) then
874  output_fields(out_num)%num_axes = 1
875  else
876  output_fields(out_num)%num_axes = 2
877  output_fields(out_num)%axes(2) = diag_axis_init('Levels',level_values,'level number', 'Y' )
878  endif
879  if(need_write_axis) then
880  lat_axis = diag_axis_init('Latitude', stations(1:num_stations)%lat,'station latitudes', 'n')
881  lon_axis = diag_axis_init('Longitude',stations(1:num_stations)%lon, 'station longitudes', 'n')
882  endif
883  need_write_axis = .false.
884 
885 ! call mpp_sync()
886 
887 end function register_station_field3d
888 
889 
890 !> \brief find_output_field returns the index of the reuqested field within the
891 !! output_fields array.
892 function find_output_field(module_name, field_name)
894  character(len=*), intent(in) :: module_name !< Module name of the requested station
895  character(len=*), intent(in) :: field_name !< Field name of the requested station
896  integer :: i
897 
898  find_output_field = -1
899  do i = 1, num_output_fields
900  if(trim(output_fields(i)%module_name) == trim(module_name) .and. &
901  lowercase(trim(output_fields(i)%output_name)) == &
902  lowercase(trim(field_name))) then
904  return
905  endif
906  end do
907 end function find_output_field
908 
909 
910 !> \brief opening_file opens a file and writes axis meta_data for all files
911 !! (only on ROOT PE, do nothing on other PEs)
912 !!
913 !! \throw FATAL, "<output_name> has axis_id = -1"
914 subroutine opening_file(file)
915 ! open file, write axis meta_data for all files (only on ROOT PE,
916 ! do nothing on other PEs)
917  integer, intent(in) :: file !< File to be opened
918  character(len=128) :: time_units
919  integer :: j,field_num,num_axes,axes(5),k
920  logical :: time_ops
921  integer :: time_axis_id(1),time_bounds_id(1)
922 
923  write(time_units, 11) trim(time_unit_list(files(file)%time_units)), base_year, &
925 11 format(a, ' since ', i4.4, '-', i2.2, '-', i2.2, ' ', i2.2, ':', i2.2, ':', i2.2)
926  call mpp_open(files(file)%file_unit, files(file)%name, action=mpp_overwr, &
927  form=mpp_netcdf, threading=mpp_single, fileset=mpp_single)
928  call mpp_write_meta (files(file)%file_unit, 'title', cval=trim(global_descriptor))
929  time_ops = .false.
930  do j = 1, files(file)%num_fields
931  field_num = files(file)%fields(j)
932  if(output_fields(field_num)%time_ops) then
933  time_ops = .true.
934  exit
935  endif
936  enddo
937 !write axis meta data
938  do j = 1, files(file)%num_fields
939  field_num = files(file)%fields(j)
940  num_axes = output_fields(field_num)%num_axes
941  axes(1:num_axes) = output_fields(field_num)%axes(1:num_axes)
942  do k = 1,num_axes
943  if(axes(k)<0) &
944  call error_mesg ('station_data opening_file','output_name '// &
945  trim(output_fields(field_num)%output_name)// &
946  ' has axis_id = -1', fatal)
947  enddo
948  axes(num_axes + 1) = lat_axis
949  axes(num_axes + 2) = lon_axis
950  axes(num_axes + 3) = files(file)%time_axis_id
951 ! need in write_axis: name, unit,long_name
952  call write_axis_meta_data(files(file)%file_unit, axes(1:num_axes + 3), time_ops)
953  if(time_ops) then
954  axes(num_axes + 4) = files(file)%time_bounds_id
955  call write_axis_meta_data(files(file)%file_unit, axes(1:num_axes + 4))
956  endif
957  end do
958 !write field meta data
959  do j = 1, files(file)%num_fields
960  field_num = files(file)%fields(j)
961  num_axes = output_fields(field_num)%num_axes
962  axes(1:num_axes) = output_fields(field_num)%axes(1:num_axes)
963  num_axes = num_axes + 1
964  axes(num_axes) = files(file)%time_axis_id
965  diag_field = write_field_meta_data(files(file)%file_unit,output_fields(field_num)%output_name, &
966  axes(1:num_axes),output_fields(field_num)%units, &
967  output_fields(field_num)%long_name,time_method=output_fields(field_num)%time_method,&
968  pack=output_fields(field_num)%pack)
969  output_fields(field_num)%f_type = diag_field%Field
970  end do
971  if(time_ops) then
972  time_axis_id(1) = files(file)%time_axis_id
973  time_bounds_id(1) = files(file)%time_bounds_id
974  diag_field=write_field_meta_data(files(file)%file_unit, avg_name // '_T1',time_axis_id, &
975  time_units,"Start time for average period", pack=1)
976  files(file)%f_avg_start = diag_field%Field
977 
978  diag_field=write_field_meta_data(files(file)%file_unit,avg_name // '_T2' ,time_axis_id, &
979  time_units,"End time for average period", pack=1)
980  files(file)%f_avg_end = diag_field%Field
981 
982  diag_field=write_field_meta_data(files(file)%file_unit,avg_name // '_DT' ,time_axis_id, &
983  time_units,"Length of average period", pack=1)
984  files(file)%f_avg_nitems = diag_field%Field
985 
986  diag_field=write_field_meta_data(files(file)%file_unit, 'Time_bnds', (/time_bounds_id,time_axis_id/), &
987  trim(time_unit_list(files(file)%time_units)), &
988  'Time axis boundaries', pack=1)
989  files(file)%f_bounds = diag_field%Field
990  endif
991  call done_meta_data(files(file)%file_unit)
992 end subroutine opening_file
993 
994 
995 !> \brief send_station_data_2d sends data to the root PE, which then sends
996 !! data to staton_data_out to be sent to files.<br>
997 !! send_station_data_2d is overloaded by send_station_data_3d as a
998 !! part of the send_station_data interface.
999 subroutine send_station_data_2d(field_id, data, time)
1000  integer, intent(in) :: field_id !< Field id of the station data being recorded
1001  real, intent(in) :: data(:,:) !< Data of the station being recorded
1002  type(time_type), intent(in) :: time !< Time of the station being recorded
1003  real :: data3d(size(data,1),size(data,2),1)
1004 
1005  data3d(:,:,1) = data
1006  call send_station_data_3d(field_id, data3d, time)
1007 end subroutine send_station_data_2d
1008 
1009 
1010 !> \brief send_station_data_3d sends data to the root PE, which then sends
1011 !! data to staton_data_out to be sent to files.<br>
1012 !! send_station_data_3d is overloaded by send_station_data_2d as a
1013 !! part of the send_station_data interface.
1014 !!
1015 !! \throw FATAL, "Station data NOT initialized"
1016 !! \throw FATAL, "counter=0 for averaged field"
1017 !! \throw FATAL, "Global field contains MISSING field"
1018 !! \throw FATAL, "Local index out of range for field"
1019 subroutine send_station_data_3d(field_id, data, time)
1021  integer, intent(in) :: field_id !< Field id of the station data being recorded
1022  real, intent(in) :: data(:,:,:) !< Data of the station being recorded
1023  type(time_type), intent(in) :: time !< Time of the station being recorded
1024  integer :: freq,units,file_num,local_num_stations,i,ii, max_counter
1025  integer :: index_x, index_y, station_id
1026  integer, allocatable :: station_ids(:) !< root_pe only, to receive local station_ids
1027  real, allocatable :: tmp_buffer(:,:) !< root_pe only, to receive local buffer from each PE
1028 
1029 
1030  if (.not.module_is_initialized) &
1031  call error_mesg ('send_station_data_3d',' station_data NOT initialized', fatal)
1032 
1033  if(field_id < 0) return
1034  file_num = output_fields(field_id)%output_file
1035  if( mpp_pe() == mpp_root_pe() .and. files(file_num)%file_unit < 0) then
1036  call opening_file(file_num)
1037  endif
1038  freq = files(file_num)%output_freq
1039  units = files(file_num)%output_units
1040  !> compare time with next_output
1041 
1042  if (time > output_fields(field_id)%next_output .and. freq /= end_of_run) then ! time to write out
1043 ! ALL PEs, including root PE, must send data to root PE
1044  call mpp_send(output_fields(field_id)%num_station,plen=1,to_pe=mpp_root_pe(),tag=comm_tag_1)
1045  if(output_fields(field_id)%num_station > 0) then
1046  call mpp_send(output_fields(field_id)%station_id(1),plen=size(output_fields(field_id)%station_id),&
1047  to_pe=mpp_root_pe())
1048  call mpp_send(output_fields(field_id)%buffer(1,1),plen=size(output_fields(field_id)%buffer),&
1049  to_pe=mpp_root_pe())
1050  endif
1051  !> get max_counter if the field is averaged
1052  if(output_fields(field_id)%time_average) then
1053  max_counter = output_fields(field_id)%counter
1054  call mpp_max(max_counter, pelist)
1055  endif
1056  !> receive local data from all PEs
1057  if(mpp_pe() == mpp_root_pe()) then
1058  do i = 1,size(pelist)
1059  call mpp_recv(local_num_stations,glen=1,from_pe=pelist(i),tag=comm_tag_1)
1060  if(local_num_stations> 0) then
1061  allocate(station_ids(local_num_stations))
1062  allocate(tmp_buffer(local_num_stations,output_fields(field_id)%nlevel))
1063  call mpp_recv(station_ids(1), glen=size(station_ids), from_pe=pelist(i))
1064  call mpp_recv(tmp_buffer(1,1),glen=size(tmp_buffer), from_pe=pelist(i))
1065  do ii = 1,local_num_stations
1066  global_field%buffer(station_ids(ii),:) = tmp_buffer(ii,:)
1067  enddo
1068  deallocate(station_ids, tmp_buffer)
1069  endif
1070  enddo
1071  !> send global_buffer content to file
1072  if(output_fields(field_id)%time_average) then
1073  if(max_counter == 0 ) &
1074  call error_mesg ('send_station_data','counter=0 for averaged field '// &
1075  output_fields(field_id)%output_name, fatal)
1076  global_field%buffer = global_field%buffer/real(max_counter)
1077  endif
1078  !> check if global_field contains any missing values
1079  if(any(global_field%buffer == missing)) &
1080  call error_mesg ('send_station_data','Global_field contains MISSING, field '// &
1081  output_fields(field_id)%output_name, fatal)
1082  call station_data_out(file_num,field_id,global_field%buffer,output_fields(field_id)%next_output)
1083  global_field%buffer = missing
1084  endif
1085  call mpp_sync()
1086  !> clear buffer, increment next_output time and reset counter on ALL PEs
1087  if(output_fields(field_id)%num_station>0) output_fields(field_id)%buffer = empty
1088  output_fields(field_id)%last_output = output_fields(field_id)%next_output
1089  output_fields(field_id)%next_output = diag_time_inc(output_fields(field_id)%next_output,&
1090  freq, units)
1091  output_fields(field_id)%counter = 0; max_counter = 0
1092 
1093  endif
1094  !> accumulate buffer only
1095  do i = 1 , output_fields(field_id)%num_station
1096  station_id = output_fields(field_id)%station_id(i)
1097  index_x = stations(station_id)%local_i; index_y = stations(station_id)%local_j
1098  if(index_x>size(data,1) .or. index_y>size(data,2)) &
1099  call error_mesg ('send_station_data','local index out of range for field '// &
1100  output_fields(field_id)%output_name, fatal)
1101  if(output_fields(field_id)%time_average) then
1102  output_fields(field_id)%buffer(i,:) = output_fields(field_id)%buffer(i,:) + &
1103  data(index_x,index_y,:) ! accumulate buffer
1104  else ! not average
1105  output_fields(field_id)%buffer(i,:) = data(index_x,index_y,:)
1106  endif
1107  enddo
1108  if(output_fields(field_id)%time_average) &
1109  output_fields(field_id)%counter = output_fields(field_id)%counter + 1
1110 end subroutine send_station_data_3d
1111 
1112 
1113 !> \brief station_data_out is responsible for sending station data to files.
1114 subroutine station_data_out(file, field, data, time,final_call_in)
1116  integer, intent(in) :: file !< Integer identifying the file to be output
1117  integer, intent(in) :: field !< Integer identifying the filed to be output
1118  real, intent(inout) :: data(:, :) !< Data to be output
1119  type(time_type), intent(in) :: time !< Time of output
1120  logical, optional, intent(in):: final_call_in !< Optional logical expression
1121  logical :: final_call
1122  integer :: i, num
1123  real :: dif, time_data(2, 1, 1), dt_time(1, 1, 1), start_dif, end_dif
1124 
1125  final_call = .false.
1126  if(present(final_call_in)) final_call = final_call_in
1127  dif = get_date_dif(time, base_time, files(file)%time_units)
1128  call mpp_write(files(file)%file_unit,output_fields(field)%f_type, data, dif)
1129  start_dif = get_date_dif(output_fields(field)%last_output, base_time,files(file)%time_units)
1130  end_dif = dif
1131  do i = 1, files(file)%num_fields
1132  num = files(file)%fields(i)
1133  if(output_fields(num)%time_ops) then
1134  if(num == field) then
1135  time_data(1, 1, 1) = start_dif
1136  call mpp_write(files(file)%file_unit, files(file)%f_avg_start, &
1137  time_data(1:1,:,:), dif)
1138  time_data(2, 1, 1) = end_dif
1139  call mpp_write(files(file)%file_unit, files(file)%f_avg_end, &
1140  time_data(2:2,:,:), dif)
1141  dt_time(1, 1, 1) = end_dif - start_dif
1142  call mpp_write(files(file)%file_unit, files(file)%f_avg_nitems, &
1143  dt_time(1:1,:,:), dif)
1144 ! Include boundary variable for CF compliance
1145  call mpp_write(files(file)%file_unit, files(file)%f_bounds, &
1146  time_data(1:2,:,:), dif)
1147  exit
1148  endif
1149  end if
1150  end do
1151  if(final_call) then
1152  if(time >= files(file)%last_flush) then
1153  call mpp_flush(files(file)%file_unit)
1154  files(file)%last_flush = time
1155  endif
1156  else
1157  if(time > files(file)%last_flush) then
1158  call mpp_flush(files(file)%file_unit)
1159  files(file)%last_flush = time
1160  endif
1161  endif
1162 end subroutine station_data_out
1163 
1164 
1165 !> \brief Must be called <tt>after the last time step</tt> to write the buffer
1166 !! content
1167 !!
1168 !! \throw FATAL, "counter=0 for averaged field"
1169 !! \throw FATAL, "Global field contains MISSING field"
1170 subroutine station_data_end(time)
1172  type(time_type), intent(in) :: time !< model's time
1173  integer :: freq, max_counter, local_num_stations
1174  integer :: file, nfield, field, pe, col
1175  integer, allocatable :: station_ids(:) !< root_pe only, to receive local station_ids
1176  real, allocatable :: tmp_buffer(:,:) !< root_pe only, to receive local buffer from each PE
1177 
1178  do file = 1, num_files
1179  freq = files(file)%output_freq
1180  do nfield = 1, files(file)%num_fields
1181  field = files(file)%fields(nfield)
1182  if(.not. output_fields(field)%register) cycle
1183  if(time >= output_fields(field)%next_output .or. freq == end_of_run) then
1184  !> ALL PEs, including root PE, must send data to root PE
1185  call mpp_send(output_fields(field)%num_station,plen=1,to_pe=mpp_root_pe(),tag=comm_tag_2)
1186  if(output_fields(field)%num_station > 0) then
1187  call mpp_send(output_fields(field)%station_id(1),plen=size(output_fields(field)%station_id),&
1188  to_pe=mpp_root_pe(),tag=comm_tag_3)
1189  call mpp_send(output_fields(field)%buffer(1,1),plen=size(output_fields(field)%buffer),&
1190  to_pe=mpp_root_pe(),tag=comm_tag_4)
1191  endif
1192  !> get max_counter if the field is averaged
1193  if(output_fields(field)%time_average) then
1194  max_counter = output_fields(field)%counter
1195  call mpp_max(max_counter, pelist)
1196  endif
1197  !> only root PE receives local data from all PEs
1198  if(mpp_pe() == mpp_root_pe()) then
1199  do pe = 1,size(pelist)
1200  call mpp_recv(local_num_stations,glen=1,from_pe=pelist(pe),tag=comm_tag_2)
1201  if(local_num_stations> 0) then
1202  allocate(station_ids(local_num_stations))
1203  allocate(tmp_buffer(local_num_stations,output_fields(field)%nlevel))
1204  call mpp_recv(station_ids(1), glen=size(station_ids), from_pe=pelist(pe),tag=comm_tag_3)
1205  call mpp_recv(tmp_buffer(1,1),glen=size(tmp_buffer),from_pe=pelist(pe),tag=comm_tag_4)
1206  do col = 1,local_num_stations
1207  global_field%buffer(station_ids(col),:) = tmp_buffer(col,:)
1208  enddo
1209  deallocate(station_ids, tmp_buffer)
1210  endif
1211  enddo
1212  !> send global_buffer content to file
1213  if(output_fields(field)%time_average)then
1214  if(max_counter == 0 )&
1215  call error_mesg ('send_station_end','counter=0 for averaged field '// &
1216  output_fields(field)%output_name, fatal)
1217  global_field%buffer = global_field%buffer/real(max_counter)
1218  endif
1219  !> check if global_field contains any missing values
1220  if(any(global_field%buffer == missing)) &
1221  call error_mesg ('send_station_end','Global_field contains MISSING, field '// &
1222  output_fields(field)%output_name, fatal)
1223  call station_data_out(file,field,global_field%buffer,output_fields(field)%next_output,.true.)
1224  global_field%buffer = missing
1225  endif
1226  call mpp_sync()
1227  endif
1228  !> deallocate field buffer
1229  if(output_fields(field)%num_station>0) &
1230  deallocate(output_fields(field)%buffer, output_fields(field)%station_id)
1231  enddo ! nfield
1232  enddo ! file
1233  if(mpp_pe() == mpp_root_pe()) deallocate(global_field%buffer)
1234 end subroutine station_data_end
1235 
1236 
1237 end module station_data_mod
1238 
1239 
Definition: fms.F90:20
type(time_type) function, public increment_time(Time, seconds, days, ticks, err_msg, allow_neg_inc)
type(time_type) function, public increment_date(Time, years, months, days, hours, minutes, seconds, ticks, err_msg, allow_neg_inc)
integer max_output_fields
subroutine, public station_data_init()
Read in lat and lon of each station Create station_id based on lat, lon Read station_data_table, initialize output_fields and output files.
integer, parameter max_files
subroutine, public done_meta_data(file_unit)
type(global_field_type), save global_field
integer function get_station_id(lat, lon)
get_station_id is passed the station&#39;s distinct lat and lon to determine what the station&#39;s id is...
integer function register_station_field3d(module_name, fieldname, glo_lat, glo_lon, levels, init_time, domain, longname, units)
register_station_field3d registers a new station field with user input. register_station_field3d is ...
integer function find_output_field(module_name, field_name)
find_output_field returns the index of the reuqested field within the output_fields array...
integer function, public nearest_index(value, array)
Definition: axis_utils.F90:443
subroutine send_station_data_2d(field_id, data, time)
send_station_data_2d sends data to the root PE, which then sends data to staton_data_out to be sent t...
integer function find_file(name)
find_file finds the index of a requested file
integer, parameter every_time
character(len=32) pelist_name
character(len=7) avg_name
subroutine opening_file(file)
opening_file opens a file and writes axis meta_data for all files (only on ROOT PE, do nothing on other PEs)
subroutine, public write_axis_meta_data(file_unit, axes, time_ops)
integer function register_station_field2d(module_name, fieldname, glo_lat, glo_lon, init_time, domain, longname, units)
register_station_field2d registers a new station field with user input. register_station_field2d is ...
Definition: mpp.F90:39
integer function, public check_nml_error(IOSTAT, NML_NAME)
Definition: fms.F90:658
subroutine check_duplicate_output_fields()
check_duplicate_output_fields takes the data pairs (output_name and output_file) and (module_name and...
character(len=input_str_length), dimension(:), allocatable, target, public input_nml_file
Definition: mpp.F90:1378
subroutine init_output_field(module_name, field_name, file_name, time_method, pack)
init_output_field initializes output_field attributes
type(station_type), dimension(:), allocatable stations
subroutine, public station_data_end(time)
Must be called after the last time step to write the buffer content.
integer function, public get_calendar_type()
subroutine send_station_data_3d(field_id, data, time)
send_station_data_3d sends data to the root PE, which then sends data to staton_data_out to be sent t...
integer function, public diag_axis_init(name, DATA, units, cart_name, long_name, direction, set_name, edges, Domain, Domain2, DomainU, aux, req, tile_count)
Definition: diag_axis.F90:157
subroutine init_file(filename, output_freq, output_units, time_units)
init_file initializes files to write station data to. Data is formatted as number of units of measure...
character(len=256) global_descriptor
integer, parameter, public no_calendar
type(time_type) function, public diag_time_inc(time, output_freq, output_units, err_msg)
Definition: diag_util.F90:1285
This module is used for outputing model results in a list of stations (not gridded arrays)...
type(group_field_type), dimension(:), allocatable, save output_fields
type(file_type), dimension(max_files), save files
type(time_type) base_time
subroutine station_data_out(file, field, data, time, final_call_in)
station_data_out is responsible for sending station data to files.
character(len=10), dimension(6) time_unit_list
integer, dimension(:), allocatable pelist
type(diag_fieldtype) function, public write_field_meta_data(file_unit, name, axes, units, long_name, range, pack, mval, avg_name, time_method, standard_name, interp_method, attributes, num_attributes, use_UGdomain)
type(diag_fieldtype), save diag_field
integer, parameter max_fields_per_file
subroutine, public error_mesg(routine, message, level)
Definition: fms.F90:529
integer, parameter end_of_run
integer num_output_fields
logical module_is_initialized