FV3 Bundle
diag_output.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 
21 #include <fms_platform.h>
22  ! <CONTACT EMAIL="seth.underwood@noaa.gov">
23  ! Seth Underwood
24  ! </CONTACT>
25 
26  ! <OVERVIEW> <TT>diag_output_mod</TT> is an integral part of
27  ! <TT>diag_manager_mod</TT>. Its function is to write axis-meta-data,
28  ! field-meta-data and field data
29  ! </OVERVIEW>
30 
31  USE mpp_io_mod, ONLY: axistype, fieldtype, mpp_io_init, mpp_open, mpp_write_meta,&
32  & mpp_write, mpp_flush, mpp_close, mpp_get_id, mpp_wronly, mpp_overwr,&
33  & mpp_netcdf, mpp_multi, mpp_single, mpp_io_unstructured_write
37  & OPERATOR(.NE.), mpp_get_layout, OPERATOR(.EQ.)
38  USE mpp_mod, ONLY: mpp_npes, mpp_pe
41  & get_domainug
44  USE fms_mod, ONLY: error_mesg, mpp_pe, write_version_number, fms_error_handler, fatal
45 
46 #ifdef use_netCDF
47  USE netcdf, ONLY: nf90_int, nf90_float, nf90_char
48 #endif
49 
50  use mpp_domains_mod, only: mpp_get_ug_io_domain
51  use mpp_domains_mod, only: mpp_get_ug_domain_npes
52  use mpp_domains_mod, only: mpp_get_ug_domain_pelist
53  use mpp_mod, only: mpp_gather
54  use mpp_mod, only: uppercase
55 
56  IMPLICIT NONE
57 
58  PRIVATE
61 
63 
64  INTEGER, PARAMETER :: netcdf1 = 1
65  INTEGER, PARAMETER :: mxch = 128
66  INTEGER, PARAMETER :: mxchl = 256
67  INTEGER :: current_file_unit = -1
68  INTEGER, DIMENSION(2,2) :: max_range = reshape((/ -32767, 32767, -127, 127 /),(/2,2/))
69 ! DATA max_range / -32767, 32767, -127, 127 /
70  INTEGER, DIMENSION(2) :: missval = (/ -32768, -128 /)
71 
72  INTEGER, PARAMETER :: max_axis_num = 20
73  INTEGER :: num_axis_in_file = 0
74  INTEGER, DIMENSION(max_axis_num) :: axis_in_file
75  LOGICAL, DIMENSION(max_axis_num) :: time_axis_flag, edge_axis_flag
76  TYPE(axistype), DIMENSION(max_axis_num), SAVE :: axis_types
77 
78  LOGICAL :: module_is_initialized = .false.
79 
80  ! Include variable "version" to be written to log file.
81 #include<file_version.h>
82 
83 CONTAINS
84 
85  ! <SUBROUTINE NAME="diag_output_init">
86  ! <OVERVIEW>
87  ! Registers the time axis and opens the output file.
88  ! </OVERVIEW>
89  ! <TEMPLATE>
90  ! SUBROUTINE diag_output_init (file_name, format, file_title, file_unit,
91  ! all_scalar_or_1d, domain)
92  ! </TEMPLATE>
93  ! <DESCRIPTION>
94  ! Registers the time axis, and opens the file for output.
95  ! </DESCRIPTION>
96  ! <IN NAME="file_name" TYPE="CHARACTER(len=*)">Output file name</IN>
97  ! <IN NAME="format" TYPE="INTEGER">File format (Currently only 'NETCDF' is valid)</IN>
98  ! <IN NAME="file_title" TYPE="CHARACTER(len=*)">Descriptive title for the file</IN>
99  ! <OUT NAME="file_unit" TYPE="INTEGER">
100  ! File unit number assigned to the output file. Needed for subsuquent calls to
101  ! <TT>diag_output_mod</TT>
102  ! </OUT>
103  ! <IN NAME="all_scalar_or_1d" TYPE="LOGICAL" />
104  ! <IN NAME="domain" TYPE="TYPE(domain2d)" />
105  ! <IN NAME="domainU" TYPE="TYPE(domainUG)" />The unstructure domain </IN>
106  SUBROUTINE diag_output_init(file_name, FORMAT, file_title, file_unit,&
107  & all_scalar_or_1d, domain, domainU, attributes)
108  CHARACTER(len=*), INTENT(in) :: file_name, file_title
109  INTEGER , INTENT(in) :: format
110  INTEGER , INTENT(out) :: file_unit
111  LOGICAL , INTENT(in) :: all_scalar_or_1d
112  TYPE(domain2d) , INTENT(in) :: domain
113  TYPE(diag_atttype), INTENT(in), DIMENSION(:), OPTIONAL :: attributes
114  TYPE(domainug), INTENT(in) :: domainu
115 
116  INTEGER :: form, threading, fileset, i
117  TYPE(diag_global_att_type) :: gatt
118 
119  !---- initialize mpp_io ----
120  IF ( .NOT.module_is_initialized ) THEN
121  CALL mpp_io_init ()
122  module_is_initialized = .true.
123  CALL write_version_number("DIAG_OUTPUT_MOD", version)
124  END IF
125 
126  !---- set up output file ----
127  SELECT CASE (format)
128  CASE (netcdf1)
129  form = mpp_netcdf
130  threading = mpp_multi
131  fileset = mpp_multi
132  CASE default
133  ! <ERROR STATUS="FATAL">invalid format</ERROR>
134  CALL error_mesg('diag_output_init', 'invalid format', fatal)
135  END SELECT
136 
137  IF(all_scalar_or_1d) THEN
138  threading = mpp_single
139  fileset = mpp_single
140  END IF
141 
142 
143 !> Check to make sure that only domain2D or domainUG is used. If both are not null, then FATAL
144  if (domain .NE. null_domain2d .AND. domainu .NE. null_domainug)&
145  & CALL error_mesg('diag_output_init', "Domain2D and DomainUG can not be used at the same time in "//&
146  & trim(file_name), fatal)
147 
148  !---- open output file (return file_unit id) -----
149  IF ( domain .NE. null_domain2d ) THEN
150  CALL mpp_open(file_unit, file_name, action=mpp_overwr, form=form,&
151  & threading=threading, fileset=fileset, domain=domain)
152  ELSE IF (domainu .NE. null_domainug) THEN
153  CALL mpp_open(file_unit, file_name, action=mpp_overwr, form=form,&
154  & threading=threading, fileset=fileset, domain_ug=domainu)
155  ELSE
156  CALL mpp_open(file_unit, file_name, action=mpp_overwr, form=form,&
157  & threading=threading, fileset=fileset)
158  END IF
159 
160  !---- write global attributes ----
161  IF ( file_title(1:1) /= ' ' ) THEN
162  CALL mpp_write_meta(file_unit, 'title', cval=trim(file_title))
163  END IF
164 
165  IF ( PRESENT(attributes) ) THEN
166  DO i=1, SIZE(attributes)
167  SELECT CASE (attributes(i)%type)
168  CASE (nf90_int)
169  CALL mpp_write_meta(file_unit, trim(attributes(i)%name), ival=attributes(i)%iatt)
170  CASE (nf90_float)
171  CALL mpp_write_meta(file_unit, trim(attributes(i)%name), rval=attributes(i)%fatt)
172  CASE (nf90_char)
173  CALL mpp_write_meta(file_unit, trim(attributes(i)%name), cval=trim(attributes(i)%catt))
174  CASE default
175  ! <ERROR STATUS="FATAL">
176  ! Unknown attribute type for attribute <name> to module/input_field <module_name>/<field_name>.
177  ! Contact the developers.
178  ! </ERROR>
179  CALL error_mesg('diag_output_mod::diag_output_init', 'Unknown attribute type for global attribute "'&
180  &//trim(attributes(i)%name)//'" in file "'//trim(file_name)//'". Contact the developers.', fatal)
181  END SELECT
182  END DO
183  END IF
184  !---- write grid type (mosaic or regular)
185  CALL get_diag_global_att(gatt)
186  CALL mpp_write_meta(file_unit, 'grid_type', cval=trim(gatt%grid_type))
187  CALL mpp_write_meta(file_unit, 'grid_tile', cval=trim(gatt%tile_name))
188 
189  END SUBROUTINE diag_output_init
190  ! </SUBROUTINE>
191 
192  ! <SUBROUTINE NAME="write_axis_meta_data">
193  ! <OVERVIEW>
194  ! Write the axes meta data to file.
195  ! </OVERVIEW>
196  ! <TEMPLATE>
197  ! SUBROUTINE write_axis_meta_data(file_unit, axes, time_ops)
198  ! </TEMPLATE>
199  ! <IN NAME="file_unit" TYPE="INTEGER">File unit number</IN>
200  ! <IN NAME="axes" TYPE="INTEGER, DIMENSION(:)">Array of axis ID's, including the time axis</IN>
201  ! <IN NAME="time_ops" TYPE="LOGICAL, OPTIONAL">
202  ! .TRUE. if this file contains any min, max, time_rms, or time_average
203  ! </IN>
204  SUBROUTINE write_axis_meta_data(file_unit, axes, time_ops)
205  INTEGER, INTENT(in) :: file_unit, axes(:)
206  LOGICAL, INTENT(in), OPTIONAL :: time_ops
207 
208  TYPE(domain1d) :: domain
209 
210  TYPE(domainug) :: domainu
211 
212  CHARACTER(len=mxch) :: axis_name, axis_units
213  CHARACTER(len=mxchl) :: axis_long_name
214  CHARACTER(len=1) :: axis_cart_name
215  INTEGER :: axis_direction, axis_edges
216  REAL, ALLOCATABLE :: axis_data(:)
217  INTEGER, ALLOCATABLE :: axis_extent(:), pelist(:)
218  INTEGER :: num_attributes
219  TYPE(diag_atttype), DIMENSION(:), ALLOCATABLE :: attributes
220  INTEGER :: calendar, id_axis, id_time_axis
221  INTEGER :: i, j, index, num, length, edges_index
222  INTEGER :: gbegin, gend, gsize, ndivs
223  LOGICAL :: time_ops1
224  CHARACTER(len=2048) :: err_msg
225  type(domainug),pointer :: io_domain
226  integer(INT_KIND) :: io_domain_npes
227  integer(INT_KIND),dimension(:),allocatable :: io_pelist
228  integer(INT_KIND),dimension(:),allocatable :: unstruct_axis_sizes
229  real,dimension(:),allocatable :: unstruct_axis_data
230 
231  ! Make sure err_msg is initialized
232  err_msg = ''
233 
234  IF ( PRESENT(time_ops) ) THEN
235  time_ops1 = time_ops
236  ELSE
237  time_ops1 = .false.
238  END IF
239 
240  !---- save the current file_unit ----
241  IF ( num_axis_in_file == 0 ) current_file_unit = file_unit
242 
243  !---- dummy checks ----
244  num = SIZE(axes(:))
245  ! <ERROR STATUS="FATAL">number of axes < 1 </ERROR>
246  IF ( num < 1 ) CALL error_mesg('write_axis_meta_data', 'number of axes < 1.', fatal)
247 
248  ! <ERROR STATUS="FATAL">writing meta data out-of-order to different files.</ERROR>
249  IF ( file_unit /= current_file_unit ) CALL error_mesg('write_axis_meta_data',&
250  & 'writing meta data out-of-order to different files.', fatal)
251 
252  !---- check all axes ----
253  !---- write axis meta data for new axes ----
254  DO i = 1, num
255  id_axis = axes(i)
256  index = get_axis_index( id_axis )
257 
258  !---- skip axes already written -----
259  IF ( index > 0 ) cycle
260 
261  !---- create new axistype (then point to) -----
263  axis_in_file(num_axis_in_file) = id_axis
265  length = get_axis_global_length(id_axis)
266  ALLOCATE(axis_data(length))
267 
268  CALL get_diag_axis(id_axis, axis_name, axis_units, axis_long_name,&
269  & axis_cart_name, axis_direction, axis_edges, domain, domainu, axis_data,&
270  & num_attributes, attributes)
271 
272  IF ( domain .NE. null_domain1d ) THEN
273  IF ( length > 0 ) THEN
274  CALL mpp_write_meta(file_unit, axis_types(num_axis_in_file),&
275  & axis_name, axis_units, axis_long_name, axis_cart_name,&
276  & axis_direction, domain, axis_data )
277  ELSE
278  CALL mpp_write_meta(file_unit, axis_types(num_axis_in_file), axis_name,&
279  & axis_units, axis_long_name, axis_cart_name, axis_direction, domain)
280  END IF
281  ELSE
282  IF ( length > 0 ) THEN
283 
284  !For an unstructured dimension, only the root rank of the io_domain
285  !pelist will perform the wirte, so a gather of the unstructured
286  !axis size and axis data is required.
287  if (uppercase(trim(axis_cart_name)) .eq. "U") then
288  if (domainu .eq. null_domainug) then
289  call error_mesg("diag_output_mod::write_axis_meta_data", &
290  "A non-nul domainUG is required to" &
291  //" write unstructured axis metadata.", &
292  fatal)
293  endif
294  io_domain => null()
295  io_domain => mpp_get_ug_io_domain(domainu)
296  io_domain_npes = mpp_get_ug_domain_npes(io_domain)
297  allocate(io_pelist(io_domain_npes))
298  call mpp_get_ug_domain_pelist(io_domain, &
299  io_pelist)
300  allocate(unstruct_axis_sizes(io_domain_npes))
301  unstruct_axis_sizes = 0
302  call mpp_gather((/size(axis_data)/), &
303  unstruct_axis_sizes, &
304  io_pelist)
305  if (mpp_pe() .eq. io_pelist(1)) then
306  allocate(unstruct_axis_data(sum(unstruct_axis_sizes)))
307  else
308  allocate(unstruct_axis_data(1))
309  endif
310  unstruct_axis_data = 0.0
311  call mpp_gather(axis_data, &
312  size(axis_data), &
313  unstruct_axis_data, &
314  unstruct_axis_sizes, &
315  io_pelist)
316  call mpp_write_meta(file_unit, &
318  axis_name, &
319  axis_units, &
320  axis_long_name, &
321  axis_cart_name, &
322  axis_direction, &
323  data=unstruct_axis_data)
324  deallocate(io_pelist)
325  deallocate(unstruct_axis_sizes)
326  deallocate(unstruct_axis_data)
327  io_domain => null()
328 
329  else
330  CALL mpp_write_meta(file_unit, axis_types(num_axis_in_file), axis_name,&
331  & axis_units, axis_long_name, axis_cart_name, axis_direction, data=axis_data)
332  endif
333 
334  ELSE
335  CALL mpp_write_meta(file_unit, axis_types(num_axis_in_file), axis_name,&
336  & axis_units, axis_long_name, axis_cart_name, axis_direction)
337  END IF
338  END IF
339 
340  ! Write axis attributes
342  CALL write_attribute_meta(file_unit, id_axis, num_attributes, attributes, err_msg)
343  IF ( len_trim(err_msg) .GT. 0 ) THEN
344  CALL error_mesg('diag_output_mod::write_axis_meta_data', trim(err_msg), fatal)
345  END IF
346 
347  !---- write additional attribute (calendar_type) for time axis ----
348  !---- NOTE: calendar attribute is compliant with CF convention
349  !---- http://www.cgd.ucar.edu/cms/eaton/netcdf/CF-current.htm#cal
350  IF ( axis_cart_name == 'T' ) THEN
352  id_time_axis = mpp_get_id(axis_types(num_axis_in_file))
353  calendar = get_calendar_type()
354  CALL mpp_write_meta(file_unit, id_time_axis, 'calendar_type', cval=trim(valid_calendar_types(calendar)))
355  CALL mpp_write_meta(file_unit, id_time_axis, 'calendar', cval=trim(valid_calendar_types(calendar)))
356  IF ( time_ops1 ) THEN
357  CALL mpp_write_meta( file_unit, id_time_axis, 'bounds', cval = trim(axis_name)//'_bnds')
358  END IF
359  ELSE
361  END IF
362 
363  DEALLOCATE(axis_data)
364 
365  ! Deallocate attributes
366  IF ( ALLOCATED(attributes) ) THEN
367  DO j=1, num_attributes
368  IF ( _allocated(attributes(j)%fatt ) ) THEN
369  DEALLOCATE(attributes(j)%fatt)
370  END IF
371  IF ( _allocated(attributes(j)%iatt ) ) THEN
372  DEALLOCATE(attributes(j)%iatt)
373  END IF
374  END DO
375  DEALLOCATE(attributes)
376  END IF
377 
378  !------------- write axis containing edge information ---------------
379 
380  ! --- this axis has no edges -----
381  IF ( axis_edges <= 0 ) cycle
382 
383  ! --- was this axis edge previously defined? ---
384  id_axis = axis_edges
385  edges_index = get_axis_index(id_axis)
386  IF ( edges_index > 0 ) cycle
387 
388  ! ---- get data for axis edges ----
389  length = get_axis_global_length( id_axis )
390  ALLOCATE(axis_data(length))
391  CALL get_diag_axis(id_axis, axis_name, axis_units, axis_long_name, axis_cart_name,&
392  & axis_direction, axis_edges, domain, domainu, axis_data, num_attributes, attributes)
393 
394  ! ---- write edges attribute to original axis ----
396  & 'edges', cval=axis_name )
397 
398  ! ---- add edges index to axis list ----
399  ! ---- assume this is not a time axis ----
401  axis_in_file(num_axis_in_file) = id_axis
404 
405  ! ---- write edges axis to file ----
406  IF ( domain .NE. null_domain1d ) THEN
407  ! assume domain decomposition is irregular and loop through all prev and next
408  ! domain pointers extracting domain extents. Assume all pes are used in
409  ! decomposition
410  CALL mpp_get_global_domain(domain, begin=gbegin, end=gend, size=gsize)
411  CALL mpp_get_layout(domain, ndivs)
412  IF ( ndivs .EQ. 1 ) THEN
413  CALL mpp_write_meta(file_unit, axis_types(num_axis_in_file), axis_name,&
414  & axis_units, axis_long_name, axis_cart_name, axis_direction, data=axis_data )
415  ELSE
416  IF ( ALLOCATED(axis_extent) ) DEALLOCATE(axis_extent)
417  ALLOCATE(axis_extent(0:ndivs-1))
418  CALL mpp_get_compute_domains(domain,size=axis_extent(0:ndivs-1))
419  gend=gend+1
420  axis_extent(ndivs-1)= axis_extent(ndivs-1)+1
421  IF ( ALLOCATED(pelist) ) DEALLOCATE(pelist)
422  ALLOCATE(pelist(0:ndivs-1))
423  CALL mpp_get_pelist(domain,pelist)
424  CALL mpp_write_meta(file_unit, axis_types(num_axis_in_file),&
425  & axis_name, axis_units, axis_long_name, axis_cart_name,&
426  & axis_direction, domain, data=axis_data)
427  END IF
428  ELSE
429  CALL mpp_write_meta(file_unit, axis_types(num_axis_in_file), axis_name, axis_units,&
430  & axis_long_name, axis_cart_name, axis_direction, data=axis_data)
431  END IF
432 
433  ! Write edge axis attributes
435  CALL write_attribute_meta(file_unit, id_axis, num_attributes, attributes, err_msg)
436  IF ( len_trim(err_msg) .GT. 0 ) THEN
437  CALL error_mesg('diag_output_mod::write_axis_meta_data', trim(err_msg), fatal)
438  END IF
439 
440  DEALLOCATE (axis_data)
441  ! Deallocate attributes
442  IF ( ALLOCATED(attributes) ) THEN
443  DO j=1, num_attributes
444  IF ( _allocated(attributes(j)%fatt ) ) THEN
445  DEALLOCATE(attributes(j)%fatt)
446  END IF
447  IF ( _allocated(attributes(j)%iatt ) ) THEN
448  DEALLOCATE(attributes(j)%iatt)
449  END IF
450  END DO
451  DEALLOCATE(attributes)
452  END IF
453  END DO
454  END SUBROUTINE write_axis_meta_data
455  ! </SUBROUTINE>
456 
457  ! <FUNCTION NAME="write_field_meta_data">
458  ! <OVERVIEW>
459  ! Write the field meta data to file.
460  ! </OVERVIEW>
461  ! <TEMPLATE>
462  ! TYPE(diag_fieldtype) FUNCTION write_field_meta_data(file_unit, name, axes, units,
463  ! long_name, rnage, pack, mval, avg_name, time_method, standard_name, interp_method)
464  ! </TEMPLATE>
465  ! <DESCRIPTION>
466  ! The meta data for the field is written to the file indicated by file_unit
467  ! </DESCRIPTION>
468  ! <IN NAME="file_unit" TYPE="INTEGER">Output file unit number</IN>
469  ! <IN NAME="name" TYPE="CHARACTER(len=*)">Field name</IN>
470  ! <IN NAME="axes" TYPE="INTEGER, DIMENSION(:)">Array of axis IDs</IN>
471  ! <IN NAME="units" TYPE="CHARACTER(len=*)">Field units</IN>
472  ! <IN NAME="long_name" TYPE="CHARACTER(len=*)">Field's long name</IN>
473  ! <IN NAME="range" TYPE="REAL, DIMENSION(2), OPTIONAL">
474  ! Valid range (min, max). If min > max, the range will be ignored
475  ! </IN>
476  ! <IN NAME="pack" TYPE="INTEGER, OPTIONAL" DEFAULT="2">
477  ! Packing flag. Only valid when range specified. Valid values:
478  ! <UL>
479  ! <LI> 1 = 64bit </LI>
480  ! <LI> 2 = 32bit </LI>
481  ! <LI> 4 = 16bit </LI>
482  ! <LI> 8 = 8bit </LI>
483  ! </UL>
484  ! </IN>
485  ! <IN NAME="mval" TYPE="REAL, OPTIONAL">Missing value, must be within valid range</IN>
486  ! <IN NAME="avg_name" TYPE="CHARACTER(len=*), OPTIONAL">
487  ! Name of variable containing time averaging info
488  ! </IN>
489  ! <IN NAME="time_method" TYPE="CHARACTER(len=*), OPTIONAL">
490  ! Name of transformation applied to the time-varying data, i.e. "avg", "min", "max"
491  ! </IN>
492  ! <IN NAME="standard_name" TYPE="CHARACTER(len=*), OPTIONAL">Standard name of field</IN>
493  ! <IN NAME="interp_method" TYPE="CHARACTER(len=*), OPTIONAL" />
494  FUNCTION write_field_meta_data ( file_unit, name, axes, units, long_name, range, pack, mval,&
495  & avg_name, time_method, standard_name, interp_method, attributes, num_attributes, &
496  & use_UGdomain) result ( Field )
497  INTEGER, INTENT(in) :: file_unit, axes(:)
498  CHARACTER(len=*), INTENT(in) :: name, units, long_name
499  REAL, OPTIONAL, INTENT(in) :: range(2), mval
500  INTEGER, OPTIONAL, INTENT(in) :: pack
501  CHARACTER(len=*), OPTIONAL, INTENT(in) :: avg_name, time_method, standard_name
502  CHARACTER(len=*), OPTIONAL, INTENT(in) :: interp_method
503  TYPE(diag_atttype), DIMENSION(:), _allocatable, OPTIONAL, INTENT(in) :: attributes
504  INTEGER, OPTIONAL, INTENT(in) :: num_attributes
505  LOGICAL, OPTIONAL, INTENT(in) :: use_ugdomain
506 
507  CHARACTER(len=256) :: standard_name2
508  CHARACTER(len=1280) :: att_str
509  TYPE(diag_fieldtype) :: field
510  LOGICAL :: coord_present
511  CHARACTER(len=40) :: aux_axes(size(axes))
512  CHARACTER(len=160) :: coord_att
513  CHARACTER(len=1024) :: err_msg
514 
515  REAL :: scale, add
516  INTEGER :: i, indexx, num, ipack, np, att_len
517  LOGICAL :: use_range
518  INTEGER :: axis_indices(size(axes))
519  logical :: use_ugdomain_local
520 
521  !---- Initialize err_msg to bank ----
522  err_msg = ''
523 
524  !---- dummy checks ----
525  coord_present = .false.
526  IF( PRESENT(standard_name) ) THEN
527  standard_name2 = standard_name
528  ELSE
529  standard_name2 = 'none'
530  END IF
531 
532  use_ugdomain_local = .false.
533  if(present(use_ugdomain)) use_ugdomain_local = use_ugdomain
534 
535  num = SIZE(axes(:))
536  ! <ERROR STATUS="FATAL">number of axes < 1</ERROR>
537  IF ( num < 1 ) CALL error_mesg ( 'write_meta_data', 'number of axes < 1', fatal)
538  ! <ERROR STATUS="FATAL">writing meta data out-of-order to different files</ERROR>
539  IF ( file_unit /= current_file_unit ) CALL error_mesg ( 'write_meta_data', &
540  & 'writing meta data out-of-order to different files', fatal)
541 
542 
543  !---- check all axes for this field ----
544  !---- set up indexing to axistypes ----
545  DO i = 1, num
546  indexx = get_axis_index(axes(i))
547  !---- point to existing axistype -----
548  IF ( indexx > 0 ) THEN
549  axis_indices(i) = indexx
550  ELSE
551  ! <ERROR STATUS="FATAL">axis data not written for field</ERROR>
552  CALL error_mesg ('write_field_meta_data',&
553  & 'axis data not written for field '//trim(name), fatal)
554  END IF
555  END DO
556 
557  ! Create coordinate attribute
558  IF ( num >= 2 .OR. (num==1 .and. use_ugdomain_local) ) THEN
559  coord_att = ' '
560  DO i = 1, num
561  aux_axes(i) = get_axis_aux(axes(i))
562  IF( trim(aux_axes(i)) /= 'none' ) THEN
563  IF(len_trim(coord_att) == 0) THEN
564  coord_att = trim(aux_axes(i))
565  ELSE
566  coord_att = trim(coord_att)// ' '//trim(aux_axes(i))
567  ENDIF
568  coord_present = .true.
569  END IF
570  END DO
571  END IF
572 
573  !--------------------- write field meta data ---------------------------
574 
575  !---- select packing? ----
576  !(packing option only valid with range option)
577  IF ( PRESENT(pack) ) THEN
578  ipack = pack
579  ELSE
580  ipack = 2
581  END IF
582 
583  !---- check range ----
584  use_range = .false.
585  add = 0.0
586  scale = 1.0
587  IF ( PRESENT(range) ) THEN
588  IF ( range(2) > range(1) ) THEN
589  use_range = .true.
590  !---- set packing parameters ----
591  IF ( ipack > 2 ) THEN
592  np = ipack/4
593  add = 0.5*(range(1)+range(2))
594  scale = (range(2)-range(1)) / real(max_range(2,np)-max_range(1,np))
595  END IF
596  END IF
597  END IF
598 
599  !---- select packing? ----
600  IF ( PRESENT(mval) ) THEN
601  field%miss = mval
602  field%miss_present = .true.
603  IF ( ipack > 2 ) THEN
604  np = ipack/4
605  field%miss_pack = REAL(missval(np))*scale+add
606  field%miss_pack_present = .true.
607  ELSE
608  field%miss_pack = mval
609  field%miss_pack_present = .false.
610  END IF
611  ELSE
612  field%miss_present = .false.
613  field%miss_pack_present = .false.
614  END IF
615 
616  !------ write meta data and return fieldtype -------
617  IF ( use_range ) THEN
618  IF ( field%miss_present ) THEN
619  CALL mpp_write_meta(file_unit, field%Field,&
620  & axis_types(axis_indices(1:num)),&
621  & name, units, long_name,&
622  & range(1), range(2),&
623  & missing=field%miss_pack,&
624  & fill=field%miss_pack,&
625  & scale=scale, add=add, pack=ipack,&
626  & time_method=time_method)
627  ELSE
628  CALL mpp_write_meta(file_unit, field%Field,&
629  & axis_types(axis_indices(1:num)),&
630  & name, units, long_name,&
631  & range(1), range(2),&
632  & missing=cmor_missing_value,&
633  & fill=cmor_missing_value,&
634  & scale=scale, add=add, pack=ipack,&
635  & time_method=time_method)
636  END IF
637  ELSE
638  IF ( field%miss_present ) THEN
639  CALL mpp_write_meta(file_unit, field%Field,&
640  & axis_types(axis_indices(1:num)),&
641  & name, units, long_name,&
642  & missing=field%miss_pack,&
643  & fill=field%miss_pack,&
644  & pack=ipack, time_method=time_method)
645  ELSE
646  CALL mpp_write_meta(file_unit, field%Field,&
647  & axis_types(axis_indices(1:num)),&
648  & name, units, long_name,&
649  & missing=cmor_missing_value,&
650  & fill=cmor_missing_value,&
651  & pack=ipack, time_method=time_method)
652  END IF
653  END IF
654 
655  !---- write user defined attributes -----
656  IF ( PRESENT(num_attributes) ) THEN
657  IF ( PRESENT(attributes) ) THEN
658  IF ( num_attributes .GT. 0 .AND. _allocated(attributes) ) THEN
659  CALL write_attribute_meta(file_unit, mpp_get_id(field%Field), num_attributes, attributes, time_method, err_msg)
660  IF ( len_trim(err_msg) .GT. 0 ) THEN
661  CALL error_mesg('diag_output_mod::write_field_meta_data',&
662  & trim(err_msg)//" Contact the developers.", fatal)
663  END IF
664  ELSE
665  ! Catch some bad cases
666  IF ( num_attributes .GT. 0 .AND. .NOT._allocated(attributes) ) THEN
667  CALL error_mesg('diag_output_mod::write_field_meta_data',&
668  & 'num_attributes > 0 but attributes is not allocated for attribute '&
669  &//trim(attributes(i)%name)//' for field '//trim(name)//'. Contact the developers.', fatal)
670  ELSE IF ( num_attributes .EQ. 0 .AND. _allocated(attributes) ) THEN
671  CALL error_mesg('diag_output_mod::write_field_meta_data',&
672  & 'num_attributes == 0 but attributes is allocated for attribute '&
673  &//trim(attributes(i)%name)//' for field '//trim(name)//'. Contact the developers.', fatal)
674  END IF
675  END IF
676  ELSE
677  ! More edge error cases
678  CALL error_mesg('diag_output_mod::write_field_meta_data',&
679  & 'num_attributes present but attributes missing for attribute '&
680  &//trim(attributes(i)%name)//' for field '//trim(name)//'. Contact the developers.', fatal)
681  END IF
682  ELSE IF ( PRESENT(attributes) ) THEN
683  CALL error_mesg('diag_output_mod::write_field_meta_data',&
684  & 'attributes present but num_attributes missing for attribute '&
685  &//trim(attributes(i)%name)//' for field '//trim(name)//'. Contact the developers.', fatal)
686  END IF
687 
688 
689  !---- write additional attribute for time averaging -----
690  IF ( PRESENT(avg_name) ) THEN
691  IF ( avg_name(1:1) /= ' ' ) THEN
692  CALL mpp_write_meta(file_unit, mpp_get_id(field%Field),&
693  & 'time_avg_info',&
694  & cval=trim(avg_name)//'_T1,'//trim(avg_name)//'_T2,'//trim(avg_name)//'_DT')
695  END IF
696  END IF
697 
698  ! write coordinates attribute for CF compliance
699  IF ( coord_present ) &
700  CALL mpp_write_meta(file_unit, mpp_get_id(field%Field),&
701  & 'coordinates', cval=trim(coord_att))
702  IF ( trim(standard_name2) /= 'none' ) CALL mpp_write_meta(file_unit, mpp_get_id(field%Field),&
703  & 'standard_name', cval=trim(standard_name2))
704 
705  !---- write attribute for interp_method ----
706  IF( PRESENT(interp_method) ) THEN
707  CALL mpp_write_meta ( file_unit, mpp_get_id(field%Field),&
708  & 'interp_method', cval=trim(interp_method))
709  END IF
710 
711  !---- get axis domain ----
712  field%Domain = get_domain2d( axes )
713  field%tile_count = get_tile_count( axes )
714  field%DomainU = get_domainug( axes(1) )
715 
716  END FUNCTION write_field_meta_data
717  ! </FUNCTION>
718 
719  !> \brief Write out attribute meta data to file
720  !!
721  !! Write out the attribute meta data to file, for field and axes
722  SUBROUTINE write_attribute_meta(file_unit, id, num_attributes, attributes, time_method, err_msg)
723  INTEGER, INTENT(in) :: file_unit !< File unit number
724  INTEGER, INTENT(in) :: id !< ID of field, file, axis to get attribute meta data
725  INTEGER, INTENT(in) :: num_attributes !< Number of attributes to write
726  TYPE(diag_atttype), DIMENSION(:), INTENT(in) :: attributes !< Array of attributes
727  CHARACTER(len=*), INTENT(in), OPTIONAL :: time_method !< To include in cell_methods attribute if present
728  CHARACTER(len=*), INTENT(out), OPTIONAL :: err_msg !< Return error message
729 
730  INTEGER :: i, att_len
731  CHARACTER(len=1280) :: att_str
732 
733  ! Clear err_msg if present
734  IF ( PRESENT(err_msg) ) err_msg = ''
735 
736  DO i = 1, num_attributes
737  SELECT CASE (attributes(i)%type)
738  CASE (nf90_int)
739  IF ( .NOT._allocated(attributes(i)%iatt) ) THEN
740  IF ( fms_error_handler('diag_output_mod::write_attribute_meta',&
741  & 'Integer attribute type indicated, but array not allocated for attribute '&
742  &//trim(attributes(i)%name)//'.', err_msg) ) THEN
743  RETURN
744  END IF
745  END IF
746  CALL mpp_write_meta(file_unit, id, trim(attributes(i)%name),&
747  & ival=attributes(i)%iatt)
748  CASE (nf90_float)
749  IF ( .NOT._allocated(attributes(i)%fatt) ) THEN
750  IF ( fms_error_handler('diag_output_mod::write_attribute_meta',&
751  & 'Real attribute type indicated, but array not allocated for attribute '&
752  &//trim(attributes(i)%name)//'.', err_msg) ) THEN
753  RETURN
754  END IF
755  END IF
756  CALL mpp_write_meta(file_unit, id, trim(attributes(i)%name),&
757  & rval=attributes(i)%fatt)
758  CASE (nf90_char)
759  att_str = attributes(i)%catt
760  att_len = attributes(i)%len
761  IF ( trim(attributes(i)%name).EQ.'cell_methods' .AND. PRESENT(time_method) ) THEN
762  ! Append ",time: time_method" if time_method present
763  att_str = attributes(i)%catt(1:attributes(i)%len)//' time: '//time_method
764  att_len = len_trim(att_str)
765  END IF
766  CALL mpp_write_meta(file_unit, id, trim(attributes(i)%name),&
767  & cval=att_str(1:att_len))
768  CASE default
769  IF ( fms_error_handler('diag_output_mod::write_attribute_meta', 'Invalid type for attribute '&
770  &//trim(attributes(i)%name)//'.', err_msg) ) THEN
771  RETURN
772  END IF
773  END SELECT
774  END DO
775  END SUBROUTINE write_attribute_meta
776 
777  ! <SUBROUTINE NAME="done_meta_data">
778  ! <OVERVIEW>
779  ! Writes axis data to file.
780  ! </OVERVIEW>
781  ! <TEMPLATE>
782  ! SUBROUTINE done_meta_data(file_unit)
783  ! </TEMPLATE>
784  ! <DESCRIPTION>
785  ! Writes axis data to file. This subroutine is to be called once per file
786  ! after all <TT>write_meta_data</TT> calls, and before the first
787  ! <TT>diag_field_out</TT> call.
788  ! </DESCRIPTION>
789  ! <IN NAME="file_unit" TYPE="INTEGER">Output file unit number</IN>
790  SUBROUTINE done_meta_data(file_unit)
791  INTEGER, INTENT(in) :: file_unit
792 
793  INTEGER :: i
794 
795  !---- write data for all non-time axes ----
796  DO i = 1, num_axis_in_file
797  IF ( time_axis_flag(i) ) cycle
798  CALL mpp_write(file_unit, axis_types(i))
799  END DO
800 
801  num_axis_in_file = 0
802  END SUBROUTINE done_meta_data
803  ! </SUBROUTINE>
804 
805  ! <SUBROUTINE NAME="diag_field_out">
806  ! <OVERVIEW>
807  ! Writes field data to an output file.
808  ! </OVERVIEW>
809  ! <TEMPLATE>
810  ! SUBROUTINE diag_field_out(file_unit, field, data, time)
811  ! </TEMPLATE>
812  ! <DESCRIPTION>
813  ! Writes field data to an output file.
814  ! </DESCRIPTION>
815  ! <IN NAME="file_unit" TYPE="INTEGER">Output file unit number</IN>
816  ! <INOUT NAME="field" TYPE="TYPE(diag_fieldtype)"></INOUT>
817  ! <INOUT NAME="data" TYPE="REAL, DIMENSIONS(:,:,:,:)"></INOUT>
818  ! <IN NAME="time" TYPE="REAL, OPTIONAL"></IN>
819  SUBROUTINE diag_field_out(file_unit, Field, DATA, time)
820  INTEGER, INTENT(in) :: file_unit
821  TYPE(diag_fieldtype), INTENT(inout) :: field
822  REAL , INTENT(inout) :: data(:,:,:,:)
823  REAL, OPTIONAL, INTENT(in) :: time
824 
825  !---- replace original missing value with (un)packed missing value ----
826  !print *, 'PE,name,miss_pack_present=',mpp_pe(), &
827  ! trim(Field%Field%name),Field%miss_pack_present
828  IF ( field%miss_pack_present ) THEN
829  WHERE ( DATA == field%miss ) DATA = field%miss_pack
830  END IF
831 
832  !---- output data ----
833  IF ( field%Domain .NE. null_domain2d ) THEN
834  IF( field%miss_present ) THEN
835  CALL mpp_write(file_unit, field%Field, field%Domain, DATA, time, &
836  tile_count=field%tile_count, default_data=field%miss_pack)
837  ELSE
838  CALL mpp_write(file_unit, field%Field, field%Domain, DATA, time, &
839  tile_count=field%tile_count, default_data=cmor_missing_value)
840  END IF
841  ELSEIF ( field%DomainU .NE. null_domainug ) THEN
842  IF( field%miss_present ) THEN
843  CALL mpp_io_unstructured_write(file_unit, field%Field, field%DomainU, DATA, tstamp=time, &
844  default_data=field%miss_pack)
845  ELSE
846  CALL mpp_io_unstructured_write(file_unit, field%Field, field%DomainU, DATA, tstamp=time, &
847  default_data=cmor_missing_value)
848  END IF
849 
850  ELSE
851  CALL mpp_write(file_unit, field%Field, DATA, time)
852  END IF
853  END SUBROUTINE diag_field_out
854  ! </SUBROUTINE>
855 
856  ! <SUBROUTINE NAME="diag_flush">
857  ! <OVERVIEW>
858  ! Flush buffer and insure data is not lost.
859  ! </OVERVIEW>
860  ! <TEMPLATE>
861  ! CALL diag_flush(file_unit)
862  ! </TEMPLATE>
863  ! <DESCRIPTION>
864  ! This subroutine can be called periodically to flush the buffer, and
865  ! insure that data is not lost if the execution fails.
866  ! </DESCRIPTION>
867  ! <IN NAME="file_unit" TYPE="INTEGER">Output file unit number to flush</IN>
868  SUBROUTINE diag_flush(file_unit)
869  INTEGER, INTENT(in) :: file_unit
870 
871  CALL mpp_flush (file_unit)
872  END SUBROUTINE diag_flush
873  ! </SUBROUTINE>
874 
875 
876  ! <FUNCTION NAME="get_axis_index">
877  ! <OVERVIEW>
878  ! Return the axis index number.
879  ! </OVERVIEW>
880  ! <TEMPLATE>
881  ! INTEGER FUNCTION get_axis_index(num)
882  ! </TEMPLATE>
883  ! <DESCRIPTION>
884  ! Return the axis index number.
885  ! </DESCRIPTION>
886  ! <IN NAME="num" TYPE="INTEGER"></IN>
887  FUNCTION get_axis_index(num) RESULT ( index )
888  INTEGER, INTENT(in) :: num
889 
890  INTEGER :: index
891  INTEGER :: i
892 
893  !---- get the array index for this axis type ----
894  !---- set up pointers to axistypes ----
895  !---- write axis meta data for new axes ----
896  index = 0
897  DO i = 1, num_axis_in_file
898  IF ( num == axis_in_file(i) ) THEN
899  index = i
900  EXIT
901  END IF
902  END DO
903  END FUNCTION get_axis_index
904  ! </FUNCTION>
905 
906  ! <SUBROUTINE NAME="get_diag_global_att">
907  ! <OVERVIEW>
908  ! Return the global attribute type.
909  ! </OVERVIEW>
910  ! <TEMPLATE>
911  ! CALL get_diag_global_att(gAtt)
912  ! </TEMPLATE>
913  ! <DESCRIPTION>
914  ! Return the global attribute type.
915  ! </DESCRIPTION>
916  ! <OUT NAME="gAtt" TYPE="TYPE(diag_global_att_type"></OUT>
917  SUBROUTINE get_diag_global_att(gAtt)
918  TYPE(diag_global_att_type), INTENT(out) :: gatt
919 
920  gatt=diag_global_att
921  END SUBROUTINE get_diag_global_att
922  ! </SUBROUTINE>
923 
924  ! <SUBROUTINE NAME="set_diag_global_att">
925  ! <OVERVIEW>
926  ! Set the global attribute type.
927  ! </OVERVIEW>
928  ! <TEMPLATE>
929  ! CALL set_diag_global_att(component, gridType, timeName)
930  ! </TEMPLATE>
931  ! <DESCRIPTION>
932  ! Set the global attribute type.
933  ! </DESCRIPTION>
934  ! <IN NAME="component" TYPE="CHARACTER(len=*)"></IN>
935  ! <IN NAME="gridType" TYPE="CHARACTER(len=*)"></IN>
936  ! <IN NAME="tileName" TYPE="CHARACTER(len=*)"></IN>
937  SUBROUTINE set_diag_global_att(component, gridType, tileName)
938  CHARACTER(len=*),INTENT(in) :: component, gridtype, tilename
939 
940  ! The following two lines are set to remove compile time warnings
941  ! about 'only used once'.
942  CHARACTER(len=64) :: component_tmp
943  component_tmp = component
944  ! Don't know how to set these for specific component
945  ! Want to be able to say
946  ! if(output_file has component) then
947  diag_global_att%grid_type = gridtype
948  diag_global_att%tile_name = tilename
949  ! endif
950  END SUBROUTINE set_diag_global_att
951  ! </SUBROUTINE>
952 
953 END MODULE diag_output_mod
Definition: fms.F90:20
real, parameter cmor_missing_value
CMOR standard missing value.
Definition: diag_data.F90:110
integer function, public get_tile_count(ids)
Definition: diag_axis.F90:799
type(diag_global_att_type), save diag_global_att
Definition: diag_output.F90:62
subroutine, public diag_field_out(file_unit, Field, DATA, time)
subroutine write_attribute_meta(file_unit, id, num_attributes, attributes, time_method, err_msg)
Write out attribute meta data to file.
type(domainug) function, public get_domainug(id)
Definition: diag_axis.F90:894
character(len=128) function, public get_axis_aux(id)
Definition: diag_axis.F90:740
type(domain1d) function, public get_domain1d(id)
Definition: diag_axis.F90:835
subroutine, public done_meta_data(file_unit)
subroutine, public diag_flush(file_unit)
integer, dimension(2, 2) max_range
Definition: diag_output.F90:68
logical function, public fms_error_handler(routine, message, err_msg)
Definition: fms.F90:573
subroutine, public diag_output_init(file_name, FORMAT, file_title, file_unit, all_scalar_or_1d, domain, domainU, attributes)
logical module_is_initialized
Definition: diag_output.F90:78
subroutine, public get_diag_axis(id, name, units, long_name, cart_name, direction, edges, Domain, DomainU, DATA, num_attributes, attributes)
Definition: diag_axis.F90:528
subroutine, public write_axis_meta_data(file_unit, axes, time_ops)
Definition: mpp.F90:39
integer function, public get_axis_global_length(id)
Definition: diag_axis.F90:778
subroutine, public set_diag_global_att(component, gridType, tileName)
integer, parameter max_axis_num
Definition: diag_output.F90:72
subroutine, public get_diag_global_att(gAtt)
integer, parameter mxch
Definition: diag_output.F90:65
integer, parameter mxchl
Definition: diag_output.F90:66
type(domain2d), save, public null_domain2d
logical, dimension(max_axis_num) edge_axis_flag
Definition: diag_output.F90:75
character(len=24) function, public valid_calendar_types(ncal, err_msg)
integer function, public get_calendar_type()
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
type(axistype), dimension(max_axis_num), save axis_types
Definition: diag_output.F90:76
integer function, public get_axis_length(id)
Definition: diag_axis.F90:712
type(domainug), save, public null_domainug
logical, dimension(max_axis_num) time_axis_flag
Definition: diag_output.F90:75
type(domain2d) function, public get_domain2d(ids)
Definition: diag_axis.F90:860
integer, dimension(max_axis_num) axis_in_file
Definition: diag_output.F90:74
integer function get_axis_index(num)
************************************************************************GNU Lesser General Public License **This file is part of the GFDL Flexible Modeling System(FMS). ! *! *FMS is free software without even the implied warranty of MERCHANTABILITY or *FITNESS FOR A PARTICULAR PURPOSE See the GNU General Public License *for more details **You should have received a copy of the GNU Lesser General Public *License along with FMS If see< http:! ***********************************************************************subroutine READ_RECORD_CORE_(unit, field, nwords, data, start, axsiz) integer, intent(in) ::unit type(fieldtype), intent(in) ::field integer, intent(in) ::nwords MPP_TYPE_, intent(inout) ::data(nwords) integer, intent(in) ::start(:), axsiz(:) integer(SHORT_KIND) ::i2vals(nwords)!rab used in conjunction with transfer intrinsic to determine size of a variable integer(KIND=1) ::one_byte(8) integer ::word_sz!#ifdef __sgi integer(INT_KIND) ::ivals(nwords) real(FLOAT_KIND) ::rvals(nwords)!#else! integer ::ivals(nwords)! real ::rvals(nwords)!#endif real(DOUBLE_KIND) ::r8vals(nwords) pointer(ptr1, i2vals) pointer(ptr2, ivals) pointer(ptr3, rvals) pointer(ptr4, r8vals) if(mpp_io_stack_size< nwords) call mpp_io_set_stack_size(nwords) call mpp_error(FATAL, 'MPP_READ currently requires use_netCDF option') end subroutine READ_RECORD_CORE_ subroutine READ_RECORD_(unit, field, nwords, data, time_level, domain, position, tile_count, start_in, axsiz_in)!routine that is finally called by all mpp_read routines to perform the read!a non-netCDF record contains:! field ID! a set of 4 coordinates(is:ie, js:je) giving the data subdomain! a timelevel and a timestamp(=NULLTIME if field is static)! 3D real data(stored as 1D)!if you are using direct access I/O, the RECL argument to OPEN must be large enough for the above!in a global direct access file, record position on PE is given by %record.!Treatment of timestamp:! We assume that static fields have been passed without a timestamp.! Here that is converted into a timestamp of NULLTIME.! For non-netCDF fields, field is treated no differently, but is written! with a timestamp of NULLTIME. There is no check in the code to prevent! the user from repeatedly writing a static field. integer, intent(in) ::unit, nwords type(fieldtype), intent(in) ::field MPP_TYPE_, intent(inout) ::data(nwords) integer, intent(in), optional ::time_level type(domain2D), intent(in), optional ::domain integer, intent(in), optional ::position, tile_count integer, intent(in), optional ::start_in(:), axsiz_in(:) integer, dimension(size(field%axes(:))) ::start, axsiz integer ::tlevel !, subdomain(4) integer ::i, error, is, ie, js, je, isg, ieg, jsg, jeg type(domain2d), pointer ::io_domain=> tlevel if(PRESENT(start_in) .AND. PRESENT(axsiz_in)) then if(size(start(! the data domain and compute domain must refer to the subdomain being passed ! In this ! since that attempts to gather all data on PE size(field%axes(:)) axsiz(i)
integer, dimension(2) missval
Definition: diag_output.F90:70
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)
integer current_file_unit
Definition: diag_output.F90:67
integer num_axis_in_file
Definition: diag_output.F90:73
subroutine, public error_mesg(routine, message, level)
Definition: fms.F90:529
integer, parameter netcdf1
Definition: diag_output.F90:64
type(domain1d), save, public null_domain1d