FMS  2024.03
Flexible Modeling System
mpp_io_read.inc
1 ! -*-f90-*-
2 
3 !***********************************************************************
4 !* GNU Lesser General Public License
5 !*
6 !* This file is part of the GFDL Flexible Modeling System (FMS).
7 !*
8 !* FMS is free software: you can redistribute it and/or modify it under
9 !* the terms of the GNU Lesser General Public License as published by
10 !* the Free Software Foundation, either version 3 of the License, or (at
11 !* your option) any later version.
12 !*
13 !* FMS is distributed in the hope that it will be useful, but WITHOUT
14 !* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
15 !* FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
16 !* for more details.
17 !*
18 !* You should have received a copy of the GNU Lesser General Public
19 !* License along with FMS. If not, see <http://www.gnu.org/licenses/>.
20 !***********************************************************************
21 !> @file
22 !> @brief Parallel file read routines for structured grids for use in @ref mpp_io_mod
23 
24 !> @addtogroup mpp_io_mod
25 !> @{
26 
27 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
28 ! !
29 ! MPP_READ !
30 ! !
31 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
32 
33 #undef MPP_READ_2DDECOMP_2D_
34 #undef READ_RECORD_CORE_
35 #define READ_RECORD_CORE_ read_record_core_r8
36 #undef READ_RECORD_
37 #define READ_RECORD_ read_record_r8
38 #define MPP_READ_2DDECOMP_2D_ mpp_read_2ddecomp_r2d_r8
39 #undef MPP_READ_2DDECOMP_3D_
40 #define MPP_READ_2DDECOMP_3D_ mpp_read_2ddecomp_r3d_r8
41 #undef MPP_READ_2DDECOMP_4D_
42 #define MPP_READ_2DDECOMP_4D_ mpp_read_2ddecomp_r4d_r8
43 #undef MPP_TYPE_
44 #define MPP_TYPE_ real(KIND=r8_kind)
45 #include <mpp_read_2Ddecomp.fh>
46 
47 #undef READ_RECORD_CORE_
48 #define READ_RECORD_CORE_ read_record_core_r4
49 #undef READ_RECORD_
50 #define READ_RECORD_ read_record_r4
51 #undef MPP_READ_2DDECOMP_2D_
52 #define MPP_READ_2DDECOMP_2D_ mpp_read_2ddecomp_r2d_r4
53 #undef MPP_READ_2DDECOMP_3D_
54 #define MPP_READ_2DDECOMP_3D_ mpp_read_2ddecomp_r3d_r4
55 #undef MPP_READ_2DDECOMP_4D_
56 #define MPP_READ_2DDECOMP_4D_ mpp_read_2ddecomp_r4d_r4
57 #undef MPP_TYPE_
58 #define MPP_TYPE_ real(KIND=r4_kind)
59 #include <mpp_read_2Ddecomp.fh>
60 
61 #undef READ_RECORD_
62 #define READ_RECORD_ read_record_r8
63 #undef MPP_READ_COMPRESSED_1D_
64 #define MPP_READ_COMPRESSED_1D_ mpp_read_compressed_r1d_r8
65 #undef MPP_READ_COMPRESSED_2D_
66 #define MPP_READ_COMPRESSED_2D_ mpp_read_compressed_r2d_r8
67 #undef MPP_READ_COMPRESSED_3D_
68 #define MPP_READ_COMPRESSED_3D_ mpp_read_compressed_r3d_r8
69 #undef MPP_TYPE_
70 #define MPP_TYPE_ real(KIND=r8_kind)
71 #include <mpp_read_compressed.fh>
72 
73 #undef READ_RECORD_
74 #define READ_RECORD_ read_record_r4
75 #undef MPP_READ_COMPRESSED_1D_
76 #define MPP_READ_COMPRESSED_1D_ mpp_read_compressed_r1d_r4
77 #undef MPP_READ_COMPRESSED_2D_
78 #define MPP_READ_COMPRESSED_2D_ mpp_read_compressed_r2d_r4
79 #undef MPP_READ_COMPRESSED_3D_
80 #define MPP_READ_COMPRESSED_3D_ mpp_read_compressed_r3d_r4
81 #undef MPP_TYPE_
82 #define MPP_TYPE_ real(KIND=r4_kind)
83 #include <mpp_read_compressed.fh>
84 
85 #include <mpp_read_distributed_ascii.inc>
86 
87 ! <SUBROUTINE NAME="mpp_read_r4D_r8" INTERFACE="mpp_read">
88 ! <IN NAME="unit" TYPE="integer"></IN>
89 ! <IN NAME="field" TYPE="type(fieldtype)"></IN>
90 ! <INOUT NAME="data" TYPE="real(KIND=r8_kind)" DIM="(:,:,:,:)"></INOUT>
91 ! <IN NAME="tindex" TYPE="integer"></IN>
92 ! </SUBROUTINE>
93  subroutine mpp_read_r4d_r8( unit, field, data, tindex)
94  integer, intent(in) :: unit
95  type(fieldtype), intent(in) :: field
96  real(KIND=r8_kind), intent(inout) :: data(:,:,:,:)
97  integer, intent(in), optional :: tindex
98 
99  call read_record_r8( unit, field, size(data(:,:,:,:)), data, tindex )
100  end subroutine mpp_read_r4d_r8
101 
102 
103 ! <SUBROUTINE NAME="mpp_read_r3D_r8" INTERFACE="mpp_read">
104 ! <IN NAME="unit" TYPE="integer"></IN>
105 ! <IN NAME="field" TYPE="type(fieldtype)"></IN>
106 ! <INOUT NAME="data" TYPE="real(KIND=r8_kind)" DIM="(:,:,:)"></INOUT>
107 ! <IN NAME="tindex" TYPE="integer"></IN>
108 ! </SUBROUTINE>
109  subroutine mpp_read_r3d_r8( unit, field, data, tindex)
110  integer, intent(in) :: unit
111  type(fieldtype), intent(in) :: field
112  real(KIND=r8_kind), intent(inout) :: data(:,:,:)
113  integer, intent(in), optional :: tindex
114 
115  call read_record_r8( unit, field, size(data(:,:,:)), data, tindex )
116  end subroutine mpp_read_r3d_r8
117 
118 ! <SUBROUTINE NAME="mpp_read_r2D_r8" INTERFACE="mpp_read">
119 ! <IN NAME="unit" TYPE="integer"></IN>
120 ! <IN NAME="field" TYPE="type(fieldtype)"></IN>
121 ! <INOUT NAME="data" TYPE="real(KIND=r8_kind)" DIM="(:,:,:)"></INOUT>
122 ! <IN NAME="tindex" TYPE="integer"></IN>
123 ! </SUBROUTINE>
124  subroutine mpp_read_r2d_r8( unit, field, data, tindex )
125  integer, intent(in) :: unit
126  type(fieldtype), intent(in) :: field
127  real(KIND=r8_kind), intent(inout) :: data(:,:)
128  integer, intent(in), optional :: tindex
129 
130  call read_record_r8( unit, field, size(data(:,:)), data, tindex )
131  end subroutine mpp_read_r2d_r8
132 
133 ! <SUBROUTINE NAME="mpp_read_r1D_r8" INTERFACE="mpp_read">
134 ! <IN NAME="unit" TYPE="integer"></IN>
135 ! <IN NAME="field" TYPE="type(fieldtype)"></IN>
136 ! <INOUT NAME="data" TYPE="real(KIND=r8_kind)" DIM="(:,:,:)"></INOUT>
137 ! <IN NAME="tindex" TYPE="integer"></IN>
138 ! </SUBROUTINE>
139  subroutine mpp_read_r1d_r8( unit, field, data, tindex )
140  integer, intent(in) :: unit
141  type(fieldtype), intent(in) :: field
142  real(KIND=r8_kind), intent(inout) :: data(:)
143  integer, intent(in), optional :: tindex
144 
145  call read_record_r8( unit, field, size(data(:)), data, tindex )
146  end subroutine mpp_read_r1d_r8
147 
148 ! <SUBROUTINE NAME="mpp_read_r0D_r8" INTERFACE="mpp_read">
149 ! <IN NAME="unit" TYPE="integer"></IN>
150 ! <IN NAME="field" TYPE="type(fieldtype)"></IN>
151 ! <INOUT NAME="data" TYPE="real(KIND=r8_kind)" DIM="(:,:,:)"></INOUT>
152 ! <IN NAME="tindex" TYPE="integer"></IN>
153 ! </SUBROUTINE>
154  subroutine mpp_read_r0d_r8( unit, field, data, tindex )
155  integer, intent(in) :: unit
156  type(fieldtype), intent(in) :: field
157  real(KIND=r8_kind), intent(inout) :: data
158  integer, intent(in), optional :: tindex
159  real(KIND=r8_kind), dimension(1) :: data_tmp
160 
161  data_tmp(1)=data
162  call read_record_r8( unit, field, 1, data_tmp, tindex )
163  data=data_tmp(1)
164  end subroutine mpp_read_r0d_r8
165 
166 ! <SUBROUTINE NAME="mpp_read_r4D_r4" INTERFACE="mpp_read">
167 ! <IN NAME="unit" TYPE="integer"></IN>
168 ! <IN NAME="field" TYPE="type(fieldtype)"></IN>
169 ! <INOUT NAME="data" TYPE="real(KIND=r8_kind)" DIM="(:,:,:,:)"></INOUT>
170 ! <IN NAME="tindex" TYPE="integer"></IN>
171 ! </SUBROUTINEr4
172  subroutine mpp_read_r4d_r4( unit, field, data, tindex)
173  integer, intent(in) :: unit
174  type(fieldtype), intent(in) :: field
175  real(KIND=r4_kind), intent(inout) :: data(:,:,:,:)
176  integer, intent(in), optional :: tindex
177 
178  call read_record_r4( unit, field, size(data(:,:,:,:)), data, tindex )
179  end subroutine mpp_read_r4d_r4
180 
181 
182 ! <SUBROUTINE NAME="mpp_read_r3D_r4" INTERFACE="mpp_read">
183 ! <IN NAME="unit" TYPE="integer"></IN>
184 ! <IN NAME="field" TYPE="type(fieldtype)"></IN>
185 ! <INOUT NAME="data" TYPE="real(KIND=r4_kind)" DIM="(:,:,:)"></INOUT>
186 ! <IN NAME="tindex" TYPE="integer"></IN>
187 ! </SUBROUTINE>
188  subroutine mpp_read_r3d_r4( unit, field, data, tindex)
189  integer, intent(in) :: unit
190  type(fieldtype), intent(in) :: field
191  real(KIND=r4_kind), intent(inout) :: data(:,:,:)
192  integer, intent(in), optional :: tindex
193 
194  call read_record_r4( unit, field, size(data(:,:,:)), data, tindex )
195  end subroutine mpp_read_r3d_r4
196 
197 ! <SUBROUTINE NAME="mpp_read_r2D_r4" INTERFACE="mpp_read">
198 ! <IN NAME="unit" TYPE="integer"></IN>
199 ! <IN NAME="field" TYPE="type(fieldtype)"></IN>
200 ! <INOUT NAME="data" TYPE="real(KIND=r4_kind)" DIM="(:,:,:)"></INOUT>
201 ! <IN NAME="tindex" TYPE="integer"></IN>
202 ! </SUBROUTINE>
203  subroutine mpp_read_r2d_r4( unit, field, data, tindex )
204  integer, intent(in) :: unit
205  type(fieldtype), intent(in) :: field
206  real(KIND=r4_kind), intent(inout) :: data(:,:)
207  integer, intent(in), optional :: tindex
208 
209  call read_record_r4( unit, field, size(data(:,:)), data, tindex )
210  end subroutine mpp_read_r2d_r4
211 
212 ! <SUBROUTINE NAME="mpp_read_r1D_r4" INTERFACE="mpp_read">
213 ! <IN NAME="unit" TYPE="integer"></IN>
214 ! <IN NAME="field" TYPE="type(fieldtype)"></IN>
215 ! <INOUT NAME="data" TYPE="real(KIND=r4_kind)" DIM="(:,:,:)"></INOUT>
216 ! <IN NAME="tindex" TYPE="integer"></IN>
217 ! </SUBROUTINE>
218  subroutine mpp_read_r1d_r4( unit, field, data, tindex )
219  integer, intent(in) :: unit
220  type(fieldtype), intent(in) :: field
221  real(KIND=r4_kind), intent(inout) :: data(:)
222  integer, intent(in), optional :: tindex
223 
224  call read_record_r4( unit, field, size(data(:)), data, tindex )
225  end subroutine mpp_read_r1d_r4
226 
227 ! <SUBROUTINE NAME="mpp_read_r0D_r4" INTERFACE="mpp_read">
228 ! <IN NAME="unit" TYPE="integer"></IN>
229 ! <IN NAME="field" TYPE="type(fieldtype)"></IN>
230 ! <INOUT NAME="data" TYPE="real(KIND=r4_kind)" DIM="(:,:,:)"></INOUT>
231 ! <IN NAME="tindex" TYPE="integer"></IN>
232 ! </SUBROUTINE>
233  subroutine mpp_read_r0d_r4( unit, field, data, tindex )
234  integer, intent(in) :: unit
235  type(fieldtype), intent(in) :: field
236  real(KIND=r4_kind), intent(inout) :: data
237  integer, intent(in), optional :: tindex
238  real(KIND=r4_kind), dimension(1) :: data_tmp
239 
240  data_tmp(1)=data
241  call read_record_r4( unit, field, 1, data_tmp, tindex )
242  data=data_tmp(1)
243  end subroutine mpp_read_r0d_r4
244 
245  subroutine mpp_read_region_r2d_r4(unit, field, data, start, nread)
246  integer, intent(in) :: unit
247  type(fieldtype), intent(in) :: field
248  real(KIND=r4_kind), intent(inout) :: data(:,:)
249  integer, intent(in) :: start(:), nread(:)
250 
251  if(size(start(:)) .NE. 4 .OR. size(nread(:)) .NE. 4) call mpp_error(fatal, &
252  "mpp_io_read.inc(mpp_read_region_r2D_r4): size of start and nread must be 4")
253 
254  if(size(data,1) .NE. nread(1) .OR. size(data,2) .NE. nread(2)) then
255  call mpp_error( fatal, 'mpp_io_read.inc(mpp_read_region_r2D_r4): size mismatch between data and nread')
256  endif
257  if(nread(3) .NE. 1 .OR. nread(4) .NE. 1) call mpp_error(fatal, &
258  "mpp_io_read.inc(mpp_read_region_r2D_r4): nread(3) and nread(4) must be 1")
259  call read_record_core_r4(unit, field, nread(1)*nread(2), data, start, nread)
260 
261  return
262 
263 
264  end subroutine mpp_read_region_r2d_r4
265 
266  subroutine mpp_read_region_r3d_r4(unit, field, data, start, nread)
267  integer, intent(in) :: unit
268  type(fieldtype), intent(in) :: field
269  real(KIND=r4_kind), intent(inout) :: data(:,:,:)
270  integer, intent(in) :: start(:), nread(:)
271 
272  if(size(start(:)) .NE. 4 .OR. size(nread(:)) .NE. 4) call mpp_error(fatal, &
273  "mpp_io_read.inc(mpp_read_region_r3D_r4): size of start and nread must be 4")
274 
275  if(size(data,1) .NE. nread(1) .OR. size(data,2) .NE. nread(2) .OR. size(data,3) .NE. nread(3) ) then
276  call mpp_error( fatal, 'mpp_io_read.inc(mpp_read_region_r3D_r4): size mismatch between data and nread')
277  endif
278  if(nread(4) .NE. 1) call mpp_error(fatal, &
279  "mpp_io_read.inc(mpp_read_region_r3D_r4): nread(4) must be 1")
280  call read_record_core_r4(unit, field, nread(1)*nread(2)*nread(3), data, start, nread)
281 
282  return
283  end subroutine mpp_read_region_r3d_r4
284 
285  subroutine mpp_read_region_r2d_r8(unit, field, data, start, nread)
286  integer, intent(in) :: unit
287  type(fieldtype), intent(in) :: field
288  real(kind=r8_kind), intent(inout) :: data(:,:)
289  integer, intent(in) :: start(:), nread(:)
290 
291  if(size(start(:)) .NE. 4 .OR. size(nread(:)) .NE. 4) call mpp_error(fatal, &
292  "mpp_io_read.inc(mpp_read_region_r2D_r8): size of start and nread must be 4")
293 
294  if(size(data,1).NE.nread(1) .OR. size(data,2).NE.nread(2)) then
295  call mpp_error( fatal, 'mpp_io_read.inc(mpp_read_region_r2D_r8): size mismatch between data and nread')
296  endif
297  if(nread(3) .NE. 1 .OR. nread(4) .NE. 1) call mpp_error(fatal, &
298  "mpp_io_read.inc(mpp_read_region_r2D_r8): nread(3) and nread(4) must be 1")
299  call read_record_core_r8(unit, field, nread(1)*nread(2), data, start, nread)
300 
301  return
302  end subroutine mpp_read_region_r2d_r8
303 
304  subroutine mpp_read_region_r3d_r8(unit, field, data, start, nread)
305  integer, intent(in) :: unit
306  type(fieldtype), intent(in) :: field
307  real(kind=r8_kind), intent(inout) :: data(:,:,:)
308  integer, intent(in) :: start(:), nread(:)
309 
310  if(size(start(:)) .NE. 4 .OR. size(nread(:)) .NE. 4) call mpp_error(fatal, &
311  "mpp_io_read.inc(mpp_read_region_r3D_r8): size of start and nread must be 4")
312 
313  if(size(data,1).NE.nread(1) .OR. size(data,2).NE.nread(2) .OR. size(data,3).NE.nread(3) ) then
314  call mpp_error( fatal, 'mpp_io_read.inc(mpp_read_region_r3D_r8): size mismatch between data and nread')
315  endif
316  if(nread(4) .NE. 1) call mpp_error(fatal, &
317  "mpp_io_read.inc(mpp_read_region_r3D_r8): nread(4) must be 1")
318  call read_record_core_r8(unit, field, nread(1)*nread(2)*nread(3), data, start, nread)
319 
320  return
321  end subroutine mpp_read_region_r3d_r8
322 
323 
324  !--- Assume the text field is at most two-dimensional
325  !--- the level is always for the first dimension
326  subroutine mpp_read_text( unit, field, data, level )
327  integer, intent(in) :: unit
328  type(fieldtype), intent(in) :: field
329  character(len=*), intent(inout) :: data
330  integer, intent(in), optional :: level
331  integer :: lev, n
332  character(len=256) :: error_msg
333  integer, dimension(size(field%axes(:))) :: start, axsiz
334  character(len=len(data)) :: text
335 
336 #ifdef use_netCDF
337  if( .NOT.module_is_initialized )call mpp_error( fatal, 'READ_RECORD: must first call mpp_io_init.' )
338  if( .NOT.mpp_file(unit)%opened )call mpp_error( fatal, 'READ_RECORD: invalid unit number.' )
339  if( mpp_file(unit)%threading.EQ.mpp_single .AND. pe.NE.mpp_root_pe() )return
340 
341  if( .NOT.mpp_file(unit)%initialized ) call mpp_error( fatal, 'MPP_READ: must first call mpp_read_meta.' )
342  lev = 1
343  if(present(level)) lev = level
344 
345  if( verbose )print '(a,2i6,2i5)', 'MPP_READ: PE, unit, %id, level =', pe, unit, mpp_file(unit)%id, lev
346 
347  if( mpp_file(unit)%format.EQ.mpp_netcdf )then
348  start = 1
349  axsiz(:) = field%size(:)
350  if(len(data) < field%size(1) ) call mpp_error(fatal, &
351  'mpp_io(mpp_read_text): the first dimension size is greater than data length')
352  select case( field%ndim)
353  case(1)
354  if(lev .NE. 1) call mpp_error(fatal,'mpp_io(mpp_read_text): level should be 1 when ndim is 1')
355  case(2)
356  if(lev<1 .OR. lev > field%size(2)) then
357  write(error_msg,'(I5,"/",I5)') lev, field%size(2)
358  call mpp_error(fatal,'mpp_io(mpp_read_text): level out of range, level/max_level='//trim(error_msg))
359  end if
360  start(2) = lev
361  axsiz(2) = 1
362  case default
363  call mpp_error( fatal, 'MPP_READ: ndim of text field should be at most 2')
364  end select
365 
366  if( verbose )print '(a,2i6,i6,12i4)', 'mpp_read_text: PE, unit, num_words, start, axsiz=', &
367  & pe, unit, len(data), start, axsiz
368 
369  select case (field%type)
370  case(nf_char)
371  if(field%ndim==1) then
372  error = nf_get_var_text(mpp_file(unit)%ncid, field%id, text)
373  else
374  error = nf_get_vara_text(mpp_file(unit)%ncid, field%id, start, axsiz, text)
375  end if
376  call netcdf_err( error, mpp_file(unit), field=field )
377  do n = 1, len_trim(text)
378  if(text(n:n) == char(0) ) exit
379  end do
380  n = n-1
381  data = text(1:n)
382  case default
383  call mpp_error( fatal, 'mpp_read_text: the field type should be NF_CHAR' )
384  end select
385  else !non-netCDF
386  call mpp_error( fatal, 'Currently dont support non-NetCDF mpp read' )
387 
388  end if
389 #else
390  call mpp_error( fatal, 'mpp_read_text currently requires use_netCDF option' )
391 #endif
392  return
393  end subroutine mpp_read_text
394 
395 ! <SUBROUTINE NAME="mpp_read_meta">
396 
397 ! <OVERVIEW>
398 ! Read metadata.
399 ! </OVERVIEW>
400 ! <DESCRIPTION>
401 ! This routine is used to read the <LINK SRC="#metadata">metadata</LINK>
402 ! describing the contents of a file. Each file can contain any number of
403 ! fields, which are functions of 0-3 space axes and 0-1 time axes. (Only
404 ! one time axis can be defined per file). The basic metadata defined <LINK
405 ! SRC="#metadata">above</LINK> for <TT>axistype</TT> and
406 ! <TT>fieldtype</TT> are stored in <TT>mpp_io_mod</TT> and
407 ! can be accessed outside of <TT>mpp_io_mod</TT> using calls to
408 ! <TT>mpp_get_info</TT>, <TT>mpp_get_atts</TT>,
409 ! <TT>mpp_get_vars</TT> and
410 ! <TT>mpp_get_times</TT>.
411 ! </DESCRIPTION>
412 ! <TEMPLATE>
413 ! call mpp_read_meta(unit)
414 ! </TEMPLATE>
415 ! <IN NAME="unit" TYPE="integer"> </IN>
416 ! <NOTE>
417 ! <TT>mpp_read_meta</TT> must be called prior to <TT>mpp_read</TT>.
418 ! </NOTE>
419 ! </SUBROUTINE>
420  subroutine mpp_read_meta(unit, read_time)
421 !
422 ! read file attributes including dimension and variable attributes
423 ! and store in filetype structure. All of the file information
424 ! with the exception of the (variable) data is stored. Attributes
425 ! are supplied to the user by get_info,get_atts,get_axes and get_fields
426 !
427 ! every PE is eligible to call mpp_read_meta
428 !
429  integer, intent(in) :: unit
430  logical, intent(in), optional :: read_time ! read_time is set to false when file action is appending or writing.
431  ! This is temporary fix for MOM6 reopen_file issue.
432  integer :: ncid,ndim,nvar_total,natt,recdim,nv,nvar,len
433  integer :: error, i, j, istat, check_exist
434  integer :: type, nvdims, nvatts, dimid
435  integer, allocatable, dimension(:) :: dimids
436  character(len=128) :: name, attname, unlimname, attval, bounds_name
437  logical :: isdim, found_bounds, get_time_info
438  integer(i8_kind) :: checksumf
439  character(len=64) :: checksum_char
440  integer :: num_checksumf, last, is, k
441 
442  integer(i2_kind), allocatable :: i2vals(:)
443  integer(i4_kind), allocatable :: ivals(:)
444  real(KIND=r4_kind), allocatable :: rvals(:)
445  real(KIND=r8_kind), allocatable :: r8vals(:)
446 
447  get_time_info = .true.
448  if(present(read_time)) get_time_info = read_time
449 
450 #ifdef use_netCDF
451 
452  if( mpp_file(unit)%format.EQ.mpp_netcdf )then
453  ncid = mpp_file(unit)%ncid
454  error = nf_inq(ncid,ndim, nvar_total,&
455  natt, recdim);call netcdf_err( error, mpp_file(unit) )
456 
457 
458  mpp_file(unit)%ndim = ndim
459  mpp_file(unit)%natt = natt
460  mpp_file(unit)%recdimid = recdim
461 !
462 ! if no recdim exists, recdimid = -1
463 ! variable id of unlimdim and length
464 !
465  if( recdim.NE.-1 )then
466  error = nf_inq_dim( ncid, recdim, unlimname, mpp_file(unit)%time_level )
467  call netcdf_err( error, mpp_file(unit) )
468  error = nf_inq_varid( ncid, unlimname, mpp_file(unit)%id )
469  call netcdf_err( error, mpp_file(unit), string='Field='//unlimname )
470  else
471  mpp_file(unit)%time_level = -1 ! set to zero so mpp_get_info returns ntime=0 if no time axis present
472  endif
473 
474  allocate(mpp_file(unit)%Att(natt))
475  allocate(dimids(ndim))
476  allocate(mpp_file(unit)%Axis(ndim))
477 
478 !
479 ! initialize fieldtype and axis type
480 !
481 
482 
483  do i=1,ndim
484  mpp_file(unit)%Axis(i) = default_axis
485  enddo
486 
487  do i=1,natt
488  mpp_file(unit)%Att(i) = default_att
489  enddo
490 
491 !
492 ! assign global attributes
493 !
494  do i=1,natt
495  error=nf_inq_attname(ncid,nf_global,i,name)
496  call netcdf_err( error, mpp_file(unit), string=' Global attribute error.' )
497  error=nf_inq_att(ncid,nf_global,trim(name),type,len);call netcdf_err( error, mpp_file(unit), &
498  string=' Attribute='//name )
499  mpp_file(unit)%Att(i)%name = name
500  mpp_file(unit)%Att(i)%len = len
501  mpp_file(unit)%Att(i)%type = type
502 !
503 ! allocate space for att data and assign
504 !
505  select case (type)
506  case (nf_char)
507  if (len.gt.max_att_length) then
508  call mpp_error(note,'GLOBAL ATT too long - not reading this metadata')
509  len=7
510  mpp_file(unit)%Att(i)%len=len
511  mpp_file(unit)%Att(i)%catt = 'unknown'
512  else
513  error=nf_get_att_text(ncid,nf_global,name,mpp_file(unit)%Att(i)%catt)
514  call netcdf_err( error, mpp_file(unit), attr=mpp_file(unit)%att(i) )
515  if (verbose.and.pe == 0) print *, 'GLOBAL ATT ',trim(name),' ',mpp_file(unit)%Att(i)%catt(1:len)
516  endif
517 !
518 ! store integers in float arrays
519 !
520  case (nf_short)
521  allocate(mpp_file(unit)%Att(i)%fatt(len), stat=istat)
522  if ( istat .ne. 0 ) then
523  write(text,'(A)') istat
524  call mpp_error(fatal, &
525  & "mpp_io_mod(mpp_read_meta): Unable to allocate mpp_file%Att%fatt, NF_SHORT case. "//&
526  & "STAT = "//trim(text))
527  end if
528  allocate(i2vals(len), stat=istat)
529  if ( istat .ne. 0 ) then
530  write(text,'(A)') istat
531  call mpp_error(fatal, &
532  & "mpp_io_mod(mpp_read_meta): Unable to allocate temporary array i2vals. STAT = "&
533  & //trim(text))
534  end if
535  error=nf_get_att_int2(ncid,nf_global,name,i2vals)
536  call netcdf_err( error, mpp_file(unit), attr=mpp_file(unit)%att(i) )
537  if( verbose .and. pe == 0 )print *, 'GLOBAL ATT ',trim(name),' ',i2vals(1:len)
538  mpp_file(unit)%Att(i)%fatt(1:len)=i2vals(1:len)
539  deallocate(i2vals)
540  case (nf_int)
541  allocate(mpp_file(unit)%Att(i)%fatt(len), stat=istat)
542  if ( istat .ne. 0 ) then
543  write(text,'(A)') istat
544  call mpp_error(fatal, &
545  & "mpp_io_mod(mpp_read_meta): Unable to allocate mpp_file%Att%fatt, NF_INT case. "//&
546  & "STAT = "//trim(text))
547  end if
548  allocate(ivals(len), stat=istat)
549  if ( istat .ne. 0 ) then
550  write(text,'(A)') istat
551  call mpp_error(fatal, &
552  & "mpp_io_mod(mpp_read_meta): Unable to allocate temporary array ivals. STAT = "&
553  & //trim(text))
554  end if
555  error=nf_get_att_int(ncid,nf_global,name,ivals)
556  call netcdf_err( error, mpp_file(unit), attr=mpp_file(unit)%att(i) )
557  if( verbose .and. pe == 0 )print *, 'GLOBAL ATT ',trim(name),' ',ivals(1:len)
558  mpp_file(unit)%Att(i)%fatt(1:len)=ivals(1:len)
559  if(lowercase(trim(name)) == 'time_axis' .and. ivals(1)==0) &
560  mpp_file(unit)%time_level = -1 ! This file is an unlimited axis restart
561  deallocate(ivals)
562  case (nf_float)
563  allocate(mpp_file(unit)%Att(i)%fatt(len), stat=istat)
564  if ( istat .ne. 0 ) then
565  write(text,'(A)') istat
566  call mpp_error(fatal, &
567  & "mpp_io_mod(mpp_read_meta): Unable to allocate mpp_file%Att%fatt, NF_FLOAT case. "//&
568  & "STAT = "//trim(text))
569  end if
570  allocate(rvals(len), stat=istat)
571  if ( istat .ne. 0 ) then
572  write(text,'(A)') istat
573  call mpp_error(fatal, &
574  & "mpp_io_mod(mpp_read_meta): Unable to allocate temporary array rvals. STAT = "&
575  & //trim(text))
576  end if
577  error=nf_get_att_real(ncid,nf_global,name,rvals)
578  call netcdf_err( error, mpp_file(unit), attr=mpp_file(unit)%att(i) )
579  mpp_file(unit)%Att(i)%fatt(1:len)=rvals(1:len)
580  if( verbose .and. pe == 0)print *, 'GLOBAL ATT ',trim(name),' ',mpp_file(unit)%Att(i)%fatt(1:len)
581  deallocate(rvals)
582  case (nf_double)
583  allocate(mpp_file(unit)%Att(i)%fatt(len), stat=istat)
584  if ( istat .ne. 0 ) then
585  write(text,'(A)') istat
586  call mpp_error(fatal, &
587  & "mpp_io_mod(mpp_read_meta): Unable to allocate mpp_file%Att%fatt, NF_DOUBLE case. "//&
588  & "STAT = "//trim(text))
589  end if
590  allocate(r8vals(len), stat=istat)
591  if ( istat .ne. 0 ) then
592  write(text,'(A)') istat
593  call mpp_error(fatal, &
594  & "mpp_io_mod(mpp_read_meta): Unable to allocate temporary array r8vals. STAT = "&
595  & //trim(text))
596  end if
597  error=nf_get_att_double(ncid,nf_global,name,r8vals)
598  call netcdf_err( error, mpp_file(unit), attr=mpp_file(unit)%att(i) )
599  mpp_file(unit)%Att(i)%fatt(1:len)=r8vals(1:len)
600  if( verbose .and. pe == 0)print *, 'GLOBAL ATT ',trim(name),' ',mpp_file(unit)%Att(i)%fatt(1:len)
601  deallocate(r8vals)
602  end select
603 
604  enddo
605 !
606 ! assign dimension name and length
607 !
608  do i=1,ndim
609  error = nf_inq_dim(ncid,i,name,len);call netcdf_err( error, mpp_file(unit) )
610  mpp_file(unit)%Axis(i)%name = name
611  mpp_file(unit)%Axis(i)%len = len
612  enddo
613 
614  nvar=0
615  do i=1, nvar_total
616  error=nf_inq_var(ncid,i,name,type,nvdims,dimids,nvatts);call netcdf_err( error, mpp_file(unit) )
617  isdim=.false.
618  do j=1,ndim
619  if( trim(lowercase(name)).EQ.trim(lowercase(mpp_file(unit)%Axis(j)%name)) )isdim=.true.
620  enddo
621  if (.not.isdim) nvar=nvar+1
622  enddo
623  mpp_file(unit)%nvar = nvar
624  allocate(mpp_file(unit)%Var(nvar))
625 
626  do i=1,nvar
627  mpp_file(unit)%Var(i) = default_field
628  enddo
629 
630 !
631 ! assign dimension info
632 !
633  do i=1, nvar_total
634  error=nf_inq_var(ncid,i,name,type,nvdims,dimids,nvatts);call netcdf_err( error, mpp_file(unit) )
635  isdim=.false.
636  do j=1,ndim
637  if( trim(lowercase(name)).EQ.trim(lowercase(mpp_file(unit)%Axis(j)%name)) )isdim=.true.
638  enddo
639 
640  if( isdim )then
641  error=nf_inq_dimid(ncid,name,dimid);call netcdf_err( error, mpp_file(unit), string=' Axis='//name )
642  mpp_file(unit)%Axis(dimid)%type = type
643  mpp_file(unit)%Axis(dimid)%did = dimid
644  mpp_file(unit)%Axis(dimid)%id = i
645  mpp_file(unit)%Axis(dimid)%natt = nvatts
646  ! get axis values
647  if( i.NE.mpp_file(unit)%id )then ! non-record dims
648  select case (type)
649  case (nf_int)
650  len=mpp_file(unit)%Axis(dimid)%len
651  allocate(mpp_file(unit)%Axis(dimid)%data(len), stat=istat)
652  if ( istat .ne. 0 ) then
653  write(text,'(A)') istat
654  call mpp_error(fatal, &
655  & "mpp_io_mod(mpp_read_meta): Unable to allocate mpp_file%Axis%data, NF_INT case. "//&
656  & "STAT = "//trim(text))
657  end if
658  allocate(ivals(len), stat=istat)
659  if ( istat .ne. 0 ) then
660  write(text,'(A)') istat
661  call mpp_error(fatal, &
662  & "mpp_io_mod(mpp_read_meta): Unable to allocate temporary array ivals. STAT = "&
663  & //trim(text))
664  end if
665  error = nf_get_var_int(ncid,i,ivals)
666  call netcdf_err( error, mpp_file(unit), mpp_file(unit)%Axis(dimid) )
667  mpp_file(unit)%Axis(dimid)%data(1:len)=ivals(1:len)
668  deallocate(ivals)
669  case (nf_float)
670  len=mpp_file(unit)%Axis(dimid)%len
671  allocate(mpp_file(unit)%Axis(dimid)%data(len), stat=istat)
672  if ( istat .ne. 0 ) then
673  write(text,'(A)') istat
674  call mpp_error(fatal, "mpp_io_mod(mpp_read_meta): Unable to allocate mpp_file%Axis%data, "//&
675  & "NF_FLOAT case. STAT = "//trim(text))
676  end if
677  allocate(rvals(len), stat=istat)
678  if ( istat .ne. 0 ) then
679  write(text,'(A)') istat
680  call mpp_error(fatal, &
681  & "mpp_io_mod(mpp_read_meta): Unable to allocate temporary array rvals. STAT = "&
682  & //trim(text))
683  end if
684  error = nf_get_var_real(ncid,i,rvals)
685  call netcdf_err( error, mpp_file(unit), mpp_file(unit)%Axis(dimid) )
686  mpp_file(unit)%Axis(dimid)%data(1:len)=rvals(1:len)
687  deallocate(rvals)
688  case (nf_double)
689  len=mpp_file(unit)%Axis(dimid)%len
690  allocate(mpp_file(unit)%Axis(dimid)%data(len), stat=istat)
691  if ( istat .ne. 0 ) then
692  write(text,'(A)') istat
693  call mpp_error(fatal, "mpp_io_mod(mpp_read_meta): Unable to allocate mpp_file%Axis%data, "//&
694  & "NF_DOUBLE case. STAT = "//trim(text))
695  end if
696  allocate(r8vals(len), stat=istat)
697  if ( istat .ne. 0 ) then
698  write(text,'(A)') istat
699  call mpp_error(fatal, &
700  & "mpp_io_mod(mpp_read_meta): Unable to allocate temporary array r8vals. STAT = "&
701  & //trim(text))
702  end if
703  error = nf_get_var_double(ncid,i,r8vals)
704  call netcdf_err( error, mpp_file(unit), mpp_file(unit)%Axis(dimid) )
705  mpp_file(unit)%Axis(dimid)%data(1:len) = r8vals(1:len)
706  deallocate(r8vals)
707  case default
708  call mpp_error( fatal, 'Invalid data type for dimension' )
709  end select
710  else if(get_time_info) then
711  len = mpp_file(unit)%time_level
712  if( len > 0 ) then
713  allocate(mpp_file(unit)%time_values(len), stat=istat)
714  if ( istat .ne. 0 ) then
715  write(text,'(A)') istat
716  call mpp_error(fatal, &
717  & "mpp_io_mod(mpp_read_meta): Unable to allocate mpp_file%time_valuse. STAT = "&
718  & //trim(text))
719  end if
720  select case (type)
721  case (nf_float)
722  allocate(rvals(len), stat=istat)
723  if ( istat .ne. 0 ) then
724  write(text,'(A)') istat
725  call mpp_error(fatal, &
726  & "mpp_io_mod(mpp_read_meta): Unable to allocate temporary array rvals. STAT = "&
727  & //trim(text))
728  end if
729  !z1l read from root pe and broadcast to other processor.
730  !In the future we will modify the code if there is performance issue for very high MPI ranks.
731  if(mpp_pe()==mpp_root_pe()) then
732  error = nf_get_var_real(ncid,i,rvals)
733  call netcdf_err( error, mpp_file(unit), mpp_file(unit)%Axis(dimid) )
734  endif
735  call mpp_broadcast(rvals, len, mpp_root_pe())
736  mpp_file(unit)%time_values(1:len) = rvals(1:len)
737  deallocate(rvals)
738  case (nf_double)
739  allocate(r8vals(len), stat=istat)
740  if ( istat .ne. 0 ) then
741  write(text,'(A)') istat
742  call mpp_error(fatal, &
743  & "mpp_io_mod(mpp_read_meta): Unable to allocate temporary array r8vals. STAT = "&
744  & //trim(text))
745  end if
746  !z1l read from root pe and broadcast to other processor.
747  !In the future we will modify the code if there is performance issue for very high MPI ranks.
748  if(mpp_pe()==mpp_root_pe()) then
749  error = nf_get_var_double(ncid,i,r8vals)
750  call netcdf_err( error, mpp_file(unit), mpp_file(unit)%Axis(dimid) )
751  endif
752  call mpp_broadcast(r8vals, len, mpp_root_pe())
753  mpp_file(unit)%time_values(1:len) = r8vals(1:len)
754  deallocate(r8vals)
755  case default
756  call mpp_error( fatal, 'Invalid data type for dimension' )
757  end select
758  endif
759  endif
760  ! assign dimension atts
761  if( nvatts.GT.0 )allocate(mpp_file(unit)%Axis(dimid)%Att(nvatts))
762 
763  do j=1,nvatts
764  mpp_file(unit)%Axis(dimid)%Att(j) = default_att
765  enddo
766 
767  do j=1,nvatts
768  error=nf_inq_attname(ncid,i,j,attname);call netcdf_err( error, mpp_file(unit) )
769  error=nf_inq_att(ncid,i,trim(attname),type,len)
770  call netcdf_err( error, mpp_file(unit), string=' Attribute='//attname )
771 
772  mpp_file(unit)%Axis(dimid)%Att(j)%name = trim(attname)
773  mpp_file(unit)%Axis(dimid)%Att(j)%type = type
774  mpp_file(unit)%Axis(dimid)%Att(j)%len = len
775 
776  select case (type)
777  case (nf_char)
778  if (len.gt.max_att_length) call mpp_error(fatal,'DIM ATT too long')
779  error=nf_get_att_text(ncid,i,trim(attname),mpp_file(unit)%Axis(dimid)%Att(j)%catt);
780  call netcdf_err( error, mpp_file(unit), attr=mpp_file(unit)%Axis(dimid)%att(j) )
781  if( verbose .and. pe == 0 ) then
782  print *, 'AXIS ',trim(mpp_file(unit)%Axis(dimid)%name),' ATT ',trim(attname)
783  print *, mpp_file(unit)%Axis(dimid)%Att(j)%catt(1:len)
784  endif
785  ! store integers in float arrays
786  ! assume dimension data not packed
787  case (nf_short)
788  allocate(mpp_file(unit)%Axis(dimid)%Att(j)%fatt(len), stat=istat)
789  if ( istat .ne. 0 ) then
790  write(text,'(A)') istat
791  call mpp_error(fatal, &
792  & "mpp_io_mod(mpp_read_meta): Unable to allocate mpp_file%Axis%Att%fatt, "//&
793  & "NF_SHORT CASE. STAT = "//trim(text))
794  end if
795  allocate(i2vals(len), stat=istat)
796  if ( istat .ne. 0 ) then
797  write(text,'(A)') istat
798  call mpp_error(fatal, "mpp_io_mod(mpp_read_meta): "// &
799  "Unable to allocate temporary array i2vals. STAT = "//trim(text))
800  end if
801  error=nf_get_att_int2(ncid,i,trim(attname),i2vals);
802  call netcdf_err( error, mpp_file(unit), attr=mpp_file(unit)%Axis(dimid)%att(j) )
803  mpp_file(unit)%Axis(dimid)%Att(j)%fatt(1:len)=i2vals(1:len)
804  if( verbose .and. pe == 0 ) &
805  print *, 'AXIS ',trim(mpp_file(unit)%Axis(dimid)%name),' ATT ',&
806  & trim(attname),' ',mpp_file(unit)%Axis(dimid)%Att(j)%fatt
807  deallocate(i2vals)
808  case (nf_int)
809  allocate(mpp_file(unit)%Axis(dimid)%Att(j)%fatt(len), stat=istat)
810  if ( istat .ne. 0 ) then
811  write(text,'(A)') istat
812  call mpp_error(fatal, "mpp_io_mod(mpp_read_meta):"// &
813  " Unable to allocate mpp_file%Axis%Att%fatt, "// &
814  "NF_INT CASE. STAT = "//trim(text))
815  end if
816  allocate(ivals(len), stat=istat)
817  if ( istat .ne. 0 ) then
818  write(text,'(A)') istat
819  call mpp_error(fatal, "mpp_io_mod(mpp_read_meta): "// &
820  "Unable to allocate temporary array ivals. STAT = "//trim(text))
821  end if
822  error=nf_get_att_int(ncid,i,trim(attname),ivals);
823  call netcdf_err( error, mpp_file(unit), attr=mpp_file(unit)%Axis(dimid)%att(j) )
824  mpp_file(unit)%Axis(dimid)%Att(j)%fatt(1:len)=ivals(1:len)
825  if( verbose .and. pe == 0 ) &
826  print *, 'AXIS ',trim(mpp_file(unit)%Axis(dimid)%name),' ATT ',trim(attname), ' ',&
827  & mpp_file(unit)%Axis(dimid)%Att(j)%fatt
828  deallocate(ivals)
829  case (nf_float)
830  allocate(mpp_file(unit)%Axis(dimid)%Att(j)%fatt(len), stat=istat)
831  if ( istat .ne. 0 ) then
832  write(text,'(A)') istat
833  call mpp_error(fatal, "mpp_io_mod(mpp_read_meta): "// &
834  "Unable to allocate mpp_file%Axis%Att%fatt, "// &
835  & "NF_FLOAT CASE. STAT = "//trim(text))
836  end if
837  allocate(rvals(len), stat=istat)
838  if ( istat .ne. 0 ) then
839  write(text,'(A)') istat
840  call mpp_error(fatal, "mpp_io_mod(mpp_read_meta): "// &
841  "Unable to allocate temporary array rvals. STAT = "&
842  & //trim(text))
843  end if
844  error=nf_get_att_real(ncid,i,trim(attname),rvals);
845  call netcdf_err( error, mpp_file(unit), attr=mpp_file(unit)%Axis(dimid)%att(j) )
846  mpp_file(unit)%Axis(dimid)%Att(j)%fatt(1:len)=rvals(1:len)
847  if( verbose .and. pe == 0 ) &
848  print *, 'AXIS ',trim(mpp_file(unit)%Axis(dimid)%name),' ATT ',trim(attname),' ', &
849  & mpp_file(unit)%Axis(dimid)%Att(j)%fatt
850  deallocate(rvals)
851  case (nf_double)
852  allocate(mpp_file(unit)%Axis(dimid)%Att(j)%fatt(len), stat=istat)
853  if ( istat .ne. 0 ) then
854  write(text,'(A)') istat
855  call mpp_error(fatal, "mpp_io_mod(mpp_read_meta): "// &
856  "Unable to allocate mpp_file%Axis%Att%fatt, "//&
857  & "NF_DOUBLE CASE. STAT = "//trim(text))
858  end if
859  allocate(r8vals(len), stat=istat)
860  if ( istat .ne. 0 ) then
861  write(text,'(A)') istat
862  call mpp_error(fatal, "mpp_io_mod(mpp_read_meta): "//&
863  "Unable to allocate temporary array r8vals. STAT = "&
864  & //trim(text))
865  end if
866  error=nf_get_att_double(ncid,i,trim(attname),r8vals);
867  call netcdf_err( error, mpp_file(unit), attr=mpp_file(unit)%Axis(dimid)%att(j) )
868  mpp_file(unit)%Axis(dimid)%Att(j)%fatt(1:len)=r8vals(1:len)
869  if( verbose .and. pe == 0 ) &
870  print *, 'AXIS ',trim(mpp_file(unit)%Axis(dimid)%name),' ATT ',trim(attname),' ',&
871  & mpp_file(unit)%Axis(dimid)%Att(j)%fatt
872  deallocate(r8vals)
873  case default
874  call mpp_error( fatal, 'Invalid data type for dimension at' )
875  end select
876  ! assign pre-defined axis attributes
877  select case(trim(attname))
878  case('long_name')
879  mpp_file(unit)%Axis(dimid)%longname=mpp_file(unit)%Axis(dimid)%Att(j)%catt(1:len)
880  case('units')
881  mpp_file(unit)%Axis(dimid)%units=mpp_file(unit)%Axis(dimid)%Att(j)%catt(1:len)
882  case('cartesian_axis')
883  mpp_file(unit)%Axis(dimid)%cartesian=mpp_file(unit)%Axis(dimid)%Att(j)%catt(1:len)
884  case('calendar')
885  mpp_file(unit)%Axis(dimid)%calendar=mpp_file(unit)%Axis(dimid)%Att(j)%catt(1:len)
886  mpp_file(unit)%Axis(dimid)%calendar = lowercase(cut0(mpp_file(unit)%Axis(dimid)%calendar))
887  if (trim(mpp_file(unit)%Axis(dimid)%calendar) == 'none') &
888  mpp_file(unit)%Axis(dimid)%calendar = 'no_calendar'
889  if (trim(mpp_file(unit)%Axis(dimid)%calendar) == 'no_leap') &
890  mpp_file(unit)%Axis(dimid)%calendar = 'noleap'
891  if (trim(mpp_file(unit)%Axis(dimid)%calendar) == '365_days') &
892  mpp_file(unit)%Axis(dimid)%calendar = '365_day'
893  if (trim(mpp_file(unit)%Axis(dimid)%calendar) == '360_days') &
894  mpp_file(unit)%Axis(dimid)%calendar = '360_day'
895  case('calendar_type')
896  mpp_file(unit)%Axis(dimid)%calendar=mpp_file(unit)%Axis(dimid)%Att(j)%catt(1:len)
897  mpp_file(unit)%Axis(dimid)%calendar = lowercase(cut0(mpp_file(unit)%Axis(dimid)%calendar))
898  if (trim(mpp_file(unit)%Axis(dimid)%calendar) == 'none') &
899  mpp_file(unit)%Axis(dimid)%calendar = 'no_calendar'
900  if (trim(mpp_file(unit)%Axis(dimid)%calendar) == 'no_leap') &
901  mpp_file(unit)%Axis(dimid)%calendar = 'noleap'
902  if (trim(mpp_file(unit)%Axis(dimid)%calendar) == '365_days') &
903  mpp_file(unit)%Axis(dimid)%calendar = '365_day'
904  if (trim(mpp_file(unit)%Axis(dimid)%calendar) == '360_days') &
905  mpp_file(unit)%Axis(dimid)%calendar = '360_day'
906  case('compress')
907  mpp_file(unit)%Axis(dimid)%compressed=mpp_file(unit)%Axis(dimid)%Att(j)%catt(1:len)
908  case('positive')
909  attval = mpp_file(unit)%Axis(dimid)%Att(j)%catt(1:len)
910  if( attval.eq.'down' )then
911  mpp_file(unit)%Axis(dimid)%sense=-1
912  else if( attval.eq.'up' )then
913  mpp_file(unit)%Axis(dimid)%sense=1
914  endif
915  end select
916 
917  enddo
918  endif
919  enddo
920 
921  ! assign axis bounds
922  do j = 1, mpp_file(unit)%ndim
923  if(.not. associated(mpp_file(unit)%Axis(j)%data)) cycle
924  len = size(mpp_file(unit)%Axis(j)%data(:))
925  allocate(mpp_file(unit)%Axis(j)%data_bounds(len+1))
926  mpp_file(unit)%Axis(j)%name_bounds = 'none'
927  bounds_name = 'none'
928  found_bounds = .false.
929  do i = 1, mpp_file(unit)%Axis(j)%natt
930  if(trim(mpp_file(unit)%Axis(j)%Att(i)%name) == 'bounds' .OR. &
931  trim(mpp_file(unit)%Axis(j)%Att(i)%name) == 'edges' ) then
932  bounds_name = mpp_file(unit)%Axis(j)%Att(i)%catt
933  found_bounds = .true.
934  exit
935  endif
936  enddo
937  !-- loop through all the fields to locate bounds_name
938  if( found_bounds ) then
939  found_bounds = .false.
940  do i = 1, mpp_file(unit)%ndim
941  if(.not. associated(mpp_file(unit)%Axis(i)%data)) cycle
942  if(trim(mpp_file(unit)%Axis(i)%name) == trim(bounds_name)) then
943  found_bounds = .true.
944  if(size(mpp_file(unit)%Axis(i)%data(:)) .NE. len+1) &
945  call mpp_error(fatal, "mpp_read_meta: improperly size bounds for field "// &
946  trim(bounds_name)//" in file "// trim(mpp_file(unit)%name) )
947  mpp_file(unit)%Axis(j)%data_bounds(:) = mpp_file(unit)%Axis(i)%data(:)
948  exit
949  endif
950  enddo
951  if( .not. found_bounds ) then
952  do i=1, nvar_total
953  error=nf_inq_var(ncid,i,name,type,nvdims,dimids,nvatts);call netcdf_err( error, mpp_file(unit) )
954  if(trim(name) == trim(bounds_name)) then
955  found_bounds = .true.
956  if(nvdims .NE. 2) &
957  call mpp_error(fatal, "mpp_read_meta: field "//trim(bounds_name)//" in file "//&
958  trim(mpp_file(unit)%name)//" must be 2-D field")
959  if(mpp_file(unit)%Axis(dimids(1))%len .NE. 2) &
960  call mpp_error(fatal, "mpp_read_meta: first dimension size of field "// &
961  trim(mpp_file(unit)%Var(i)%name)//" from file "//trim(mpp_file(unit)%name)// &
962  " must be 2")
963  if(mpp_file(unit)%Axis(dimids(2))%len .NE. len) &
964  call mpp_error(fatal, "mpp_read_meta: second dimension size of field "// &
965  trim(mpp_file(unit)%Var(i)%name)//" from file "//trim(mpp_file(unit)%name)// &
966  " is not correct")
967  select case (type)
968  case (nf_int)
969  allocate(ivals(2*len), stat=istat)
970  if ( istat .ne. 0 ) then
971  write(text,'(A)') istat
972  call mpp_error(fatal, "mpp_io_mod(mpp_read_meta): Unable to allocate array ivals."//&
973  " STAT = "//trim(text))
974  end if
975  error = nf_get_var_int(ncid,i,ivals)
976  call netcdf_err( error, mpp_file(unit), string=" Field="//trim(bounds_name) )
977  mpp_file(unit)%Axis(j)%data_bounds(1:len) =ivals(1:(2*len-1):2)
978  mpp_file(unit)%Axis(j)%data_bounds(len+1) = ivals(2*len)
979  deallocate(ivals)
980  case (nf_float)
981  allocate(rvals(2*len), stat=istat)
982  if ( istat .ne. 0 ) then
983  write(text,'(A)') istat
984  call mpp_error(fatal, "mpp_io_mod(mpp_read_meta): Unable to allocate array rvals. "// &
985  " STAT = "//trim(text))
986  end if
987  error = nf_get_var_real(ncid,i,rvals)
988  call netcdf_err( error, mpp_file(unit), string=" Field="//trim(bounds_name) )
989  mpp_file(unit)%Axis(j)%data_bounds(1:len) =rvals(1:(2*len-1):2)
990  mpp_file(unit)%Axis(j)%data_bounds(len+1) = rvals(2*len)
991  deallocate(rvals)
992  case (nf_double)
993  allocate(r8vals(2*len), stat=istat)
994  if ( istat .ne. 0 ) then
995  write(text,'(A)') istat
996  call mpp_error(fatal, "mpp_io_mod(mpp_read_meta): Unable to allocate array r8vals. "//&
997  " STAT = "//trim(text))
998  end if
999  error = nf_get_var_double(ncid,i,r8vals)
1000  call netcdf_err( error, mpp_file(unit), string=" Field="//trim(bounds_name) )
1001  mpp_file(unit)%Axis(j)%data_bounds(1:len) =r8vals(1:(2*len-1):2)
1002  mpp_file(unit)%Axis(j)%data_bounds(len+1) = r8vals(2*len)
1003  deallocate(r8vals)
1004  case default
1005  call mpp_error( fatal, 'mpp_io_mod(mpp_read_meta): Invalid data type for dimension' )
1006  end select
1007  exit
1008  endif
1009  enddo
1010  endif
1011  endif
1012  if (found_bounds) then
1013  mpp_file(unit)%Axis(j)%name_bounds = trim(bounds_name)
1014  else
1015  deallocate(mpp_file(unit)%Axis(j)%data_bounds)
1016  mpp_file(unit)%Axis(j)%data_bounds =>null()
1017  endif
1018  enddo
1019 
1020 ! assign variable info
1021  nv = 0
1022  do i=1, nvar_total
1023  error=nf_inq_var(ncid,i,name,type,nvdims,dimids,nvatts);call netcdf_err( error, mpp_file(unit) )
1024 !
1025 ! is this a dimension variable?
1026 !
1027  isdim=.false.
1028  do j=1,ndim
1029  if( trim(lowercase(name)).EQ.trim(lowercase(mpp_file(unit)%Axis(j)%name)) )isdim=.true.
1030  enddo
1031 
1032  if( .not.isdim )then
1033 ! for non-dimension variables
1034  nv=nv+1
1035  if( nv.GT.mpp_file(unit)%nvar )call mpp_error(fatal,'variable index exceeds number of defined variables')
1036  mpp_file(unit)%Var(nv)%type = type
1037  mpp_file(unit)%Var(nv)%id = i
1038  mpp_file(unit)%Var(nv)%name = name
1039  mpp_file(unit)%Var(nv)%natt = nvatts
1040 ! determine packing attribute based on NetCDF variable type
1041  select case (type)
1042  case(nf_short)
1043  mpp_file(unit)%Var(nv)%pack = 4
1044  case(nf_float)
1045  mpp_file(unit)%Var(nv)%pack = 2
1046  case(nf_double)
1047  mpp_file(unit)%Var(nv)%pack = 1
1048  case (nf_int)
1049  mpp_file(unit)%Var(nv)%pack = 2
1050  case (nf_char)
1051  mpp_file(unit)%Var(nv)%pack = 1
1052  case default
1053  call mpp_error( fatal, 'Invalid variable type in NetCDF file' )
1054  end select
1055 ! assign dimension ids
1056  mpp_file(unit)%Var(nv)%ndim = nvdims
1057  allocate(mpp_file(unit)%Var(nv)%axes(nvdims))
1058  do j=1,nvdims
1059  mpp_file(unit)%Var(nv)%axes(j) = mpp_file(unit)%Axis(dimids(j))
1060  enddo
1061  allocate(mpp_file(unit)%Var(nv)%size(nvdims))
1062 
1063  do j=1,nvdims
1064  if(dimids(j).eq.mpp_file(unit)%recdimid .and. mpp_file(unit)%time_level/=-1)then
1065  mpp_file(unit)%Var(nv)%time_axis_index = j !dimids(j). z1l: Should be j.
1066  !This will cause problem when
1067  !appending to existed file.
1068  mpp_file(unit)%Var(nv)%size(j)=1 ! dimid length set to 1 here for consistency w/ mpp_write
1069  else
1070  mpp_file(unit)%Var(nv)%size(j)=mpp_file(unit)%Axis(dimids(j))%len
1071  endif
1072  enddo
1073 ! assign variable atts
1074  if( nvatts.GT.0 )allocate(mpp_file(unit)%Var(nv)%Att(nvatts))
1075 
1076  do j=1,nvatts
1077  mpp_file(unit)%Var(nv)%Att(j) = default_att
1078  enddo
1079 
1080  do j=1,nvatts
1081  error=nf_inq_attname(ncid,i,j,attname)
1082  call netcdf_err( error, mpp_file(unit), field=mpp_file(unit)%Var(nv) )
1083  error=nf_inq_att(ncid,i,attname,type,len)
1084  call netcdf_err( error, mpp_file(unit),field= mpp_file(unit)%Var(nv), string=' Attribute='//attname )
1085  mpp_file(unit)%Var(nv)%Att(j)%name = trim(attname)
1086  mpp_file(unit)%Var(nv)%Att(j)%type = type
1087  mpp_file(unit)%Var(nv)%Att(j)%len = len
1088 
1089  select case (type)
1090  case (nf_char)
1091  if (len.gt.512) call mpp_error(fatal,'VAR ATT too long')
1092  error=nf_get_att_text(ncid,i,trim(attname),mpp_file(unit)%Var(nv)%Att(j)%catt(1:len))
1093  call netcdf_err( error, mpp_file(unit), field=mpp_file(unit)%var(nv), &
1094  & attr=mpp_file(unit)%var(nv)%att(j) )
1095  if (verbose .and. pe == 0 )&
1096  print *, 'Var ',nv,' ATT ',trim(attname),' ',mpp_file(unit)%Var(nv)%Att(j)%catt(1:len)
1097 ! store integers as float internally
1098  case (nf_short)
1099  allocate(mpp_file(unit)%Var(nv)%Att(j)%fatt(len), stat=istat)
1100  if ( istat .ne. 0 ) then
1101  write(text,'(A)') istat
1102  call mpp_error(fatal, &
1103  & "mpp_io_mod(mpp_read_meta): Unable to allocate mpp_file%Var%Att%fatt, "//&
1104  & "NF_SHORT CASE. STAT = "//trim(text))
1105  end if
1106  allocate(i2vals(len), stat=istat)
1107  if ( istat .ne. 0 ) then
1108  write(text,'(A)') istat
1109  call mpp_error(fatal, &
1110  & "mpp_io_mod(mpp_read_meta): Unable to allocate temporary array i2vals. STAT = "&
1111  & //trim(text))
1112  end if
1113  error=nf_get_att_int2(ncid,i,trim(attname),i2vals)
1114  call netcdf_err( error, mpp_file(unit), field=mpp_file(unit)%var(nv), &
1115  & attr=mpp_file(unit)%var(nv)%att(j) )
1116  mpp_file(unit)%Var(nv)%Att(j)%fatt(1:len)= i2vals(1:len)
1117  if( verbose .and. pe == 0 )&
1118  print *, 'Var ',nv,' ATT ',trim(attname),' ',mpp_file(unit)%Var(nv)%Att(j)%fatt
1119  deallocate(i2vals)
1120  case (nf_int)
1121  allocate(mpp_file(unit)%Var(nv)%Att(j)%fatt(len), stat=istat)
1122  if ( istat .ne. 0 ) then
1123  write(text,'(A)') istat
1124  call mpp_error(fatal, &
1125  & "mpp_io_mod(mpp_read_meta): Unable to allocate mpp_file%Var%Att%fatt, "//&
1126  & "NF_INT CASE. STAT = "//trim(text))
1127  end if
1128  allocate(ivals(len), stat=istat)
1129  if ( istat .ne. 0 ) then
1130  write(text,'(A)') istat
1131  call mpp_error(fatal, "mpp_io_mod(mpp_read_meta): "// &
1132  "Unable to allocate temporary array ivals. STAT = "&
1133  & //trim(text))
1134  end if
1135  error=nf_get_att_int(ncid,i,trim(attname),ivals)
1136  call netcdf_err( error, mpp_file(unit), field=mpp_file(unit)%var(nv), &
1137  attr=mpp_file(unit)%var(nv)%att(j) )
1138  mpp_file(unit)%Var(nv)%Att(j)%fatt(1:len)=ivals(1:len)
1139  if( verbose .and. pe == 0 )&
1140  print *, 'Var ',nv,' ATT ',trim(attname),' ',mpp_file(unit)%Var(nv)%Att(j)%fatt
1141  deallocate(ivals)
1142  case (nf_float)
1143  allocate(mpp_file(unit)%Var(nv)%Att(j)%fatt(len), stat=istat)
1144  if ( istat .ne. 0 ) then
1145  write(text,'(A)') istat
1146  call mpp_error(fatal, "mpp_io_mod(mpp_read_meta): "// &
1147  "Unable to allocate mpp_file%Var%Att%fatt, "//&
1148  & "NF_FLOAT CASE. STAT = "//trim(text))
1149  end if
1150  allocate(rvals(len), stat=istat)
1151  if ( istat .ne. 0 ) then
1152  write(text,'(A)') istat
1153  call mpp_error(fatal, "mpp_io_mod(mpp_read_meta): "// &
1154  "Unable to allocate temporary array rvals. STAT = "&
1155  & //trim(text))
1156  end if
1157  error=nf_get_att_real(ncid,i,trim(attname),rvals)
1158  call netcdf_err( error, mpp_file(unit), field=mpp_file(unit)%var(nv), &
1159  attr=mpp_file(unit)%var(nv)%att(j) )
1160  mpp_file(unit)%Var(nv)%Att(j)%fatt(1:len)=rvals(1:len)
1161  if( verbose .and. pe == 0 )&
1162  print *, 'Var ',nv,' ATT ',trim(attname),' ',mpp_file(unit)%Var(nv)%Att(j)%fatt
1163  deallocate(rvals)
1164  case (nf_double)
1165  allocate(mpp_file(unit)%Var(nv)%Att(j)%fatt(len), stat=istat)
1166  if ( istat .ne. 0 ) then
1167  write(text,'(A)') istat
1168  call mpp_error(fatal,"mpp_io_mod(mpp_read_meta): Unable to allocate mpp_file%Var%Att%fatt, "//&
1169  & "NF_DOUBLE CASE. STAT = "//trim(text))
1170  end if
1171  allocate(r8vals(len), stat=istat)
1172  if ( istat .ne. 0 ) then
1173  write(text,'(A)') istat
1174  call mpp_error(fatal, "mpp_io_mod(mpp_read_meta): "// &
1175  "Unable to allocate temporary array r8vals. STAT = "&
1176  & //trim(text))
1177  end if
1178  error=nf_get_att_double(ncid,i,trim(attname),r8vals)
1179  call netcdf_err( error, mpp_file(unit), field=mpp_file(unit)%var(nv), &
1180  attr=mpp_file(unit)%var(nv)%att(j) )
1181  mpp_file(unit)%Var(nv)%Att(j)%fatt(1:len)=r8vals(1:len)
1182  if( verbose .and. pe == 0 ) &
1183  print *, 'Var ',nv,' ATT ',trim(attname),' ',mpp_file(unit)%Var(nv)%Att(j)%fatt
1184  deallocate(r8vals)
1185  case default
1186  call mpp_error( fatal, 'Invalid data type for variable att' )
1187  end select
1188 ! assign pre-defined field attributes
1189  select case (trim(attname))
1190  case ('long_name')
1191  mpp_file(unit)%Var(nv)%longname=mpp_file(unit)%Var(nv)%Att(j)%catt(1:len)
1192  case('units')
1193  mpp_file(unit)%Var(nv)%units=mpp_file(unit)%Var(nv)%Att(j)%catt(1:len)
1194  case('scale_factor')
1195  mpp_file(unit)%Var(nv)%scale=mpp_file(unit)%Var(nv)%Att(j)%fatt(1)
1196  case('missing')
1197  mpp_file(unit)%Var(nv)%missing=mpp_file(unit)%Var(nv)%Att(j)%fatt(1)
1198  case('missing_value')
1199  mpp_file(unit)%Var(nv)%missing=mpp_file(unit)%Var(nv)%Att(j)%fatt(1)
1200  case('_FillValue')
1201  mpp_file(unit)%Var(nv)%fill=mpp_file(unit)%Var(nv)%Att(j)%fatt(1)
1202  case('add_offset')
1203  mpp_file(unit)%Var(nv)%add=mpp_file(unit)%Var(nv)%Att(j)%fatt(1)
1204  case('packing')
1205  mpp_file(unit)%Var(nv)%pack=mpp_file(unit)%Var(nv)%Att(j)%fatt(1)
1206  case('valid_range')
1207  mpp_file(unit)%Var(nv)%min=mpp_file(unit)%Var(nv)%Att(j)%fatt(1)
1208  mpp_file(unit)%Var(nv)%max=mpp_file(unit)%Var(nv)%Att(j)%fatt(2)
1209  case('checksum')
1210  checksum_char = mpp_file(unit)%Var(nv)%Att(j)%catt
1211 ! Scan checksum attribute for , delimiter. If found implies multiple time levels.
1212  checksumf = 0
1213  num_checksumf = 1
1214 
1215  last = len_trim(checksum_char)
1216  is = index(trim(checksum_char),",") ! A value of 0 implies only 1 checksum value
1217  do while ((is > 0) .and. (is < (last-15)))
1218  is = is + scan(checksum_char(is:last), "," ) ! move starting pointer after ","
1219  num_checksumf = num_checksumf + 1
1220  enddo
1221  is =1
1222  do k = 1, num_checksumf
1223  read (checksum_char(is:is+15),'(Z16)') checksumf
1224  mpp_file(unit)%Var(nv)%checksum(k) = checksumf
1225  is = is + 17 ! Move index past the ,
1226  enddo
1227  end select
1228  enddo
1229  endif
1230  enddo ! end variable loop
1231  else
1232  call mpp_error( fatal, 'MPP READ CURRENTLY DOES NOT SUPPORT NON-NETCDF' )
1233  endif
1234 
1235  mpp_file(unit)%initialized = .true.
1236 #else
1237  call mpp_error( fatal, 'MPP_READ currently requires use_netCDF option' )
1238 #endif
1239  return
1240  end subroutine mpp_read_meta
1241 
1242 
1243  function cut0(string)
1244  character(len=256) :: cut0
1245  character(len=*), intent(in) :: string
1246  integer :: i
1247 
1248  cut0 = string
1249  i = index(string,achar(0))
1250  if(i > 0) cut0(i:i) = ' '
1251 
1252  return
1253  end function cut0
1254 
1255 
1256  subroutine mpp_get_tavg_info(unit, field, fields, tstamp, tstart, tend, tavg)
1257  implicit none
1258  integer, intent(in) :: unit
1259  type(fieldtype), intent(in) :: field
1260  type(fieldtype), intent(in), dimension(:) :: fields
1261  real, intent(inout), dimension(:) :: tstamp, tstart, tend, tavg
1262 !balaji: added because mpp_read can only read default reals
1263 ! when running with -r4 this will read a default real and then cast double
1264  real :: t_default_real
1265 
1266 
1267  integer :: n, m
1268  logical :: tavg_info_exists
1269 
1270  tavg = -1.0
1271 
1272 
1273  if (size(tstamp,1) /= size(tstart,1)) call mpp_error(fatal,&
1274  'size mismatch in mpp_get_tavg_info')
1275 
1276  if ((size(tstart,1) /= size(tend,1)) .OR. (size(tstart,1) /= size(tavg,1))) then
1277  call mpp_error(fatal,'size mismatch in mpp_get_tavg_info')
1278  endif
1279 
1280  tstart = tstamp
1281  tend = tstamp
1282 
1283  tavg_info_exists = .false.
1284 
1285 #ifdef use_netCDF
1286  do n= 1, field%natt
1287  if (field%Att(n)%type .EQ. nf_char) then
1288  if (field%Att(n)%name(1:13) == 'time_avg_info') then
1289  tavg_info_exists = .true.
1290  exit
1291  endif
1292  endif
1293  enddo
1294 #endif
1295  if (tavg_info_exists) then
1296  do n = 1, size(fields(:))
1297  if (trim(fields(n)%name) == 'average_T1') then
1298  do m = 1, size(tstart(:))
1299  call mpp_read(unit, fields(n),t_default_real, m)
1300  tstart(m) = t_default_real
1301  enddo
1302  endif
1303  if (trim(fields(n)%name) == 'average_T2') then
1304  do m = 1, size(tend(:))
1305  call mpp_read(unit, fields(n),t_default_real, m)
1306  tend(m) = t_default_real
1307  enddo
1308  endif
1309  if (trim(fields(n)%name) == 'average_DT') then
1310  do m = 1, size(tavg(:))
1311  call mpp_read(unit, fields(n),t_default_real, m)
1312  tavg(m) = t_default_real
1313  enddo
1314  endif
1315  enddo
1316 
1317  end if
1318  return
1319  end subroutine mpp_get_tavg_info
1320 !> @}
integer function mpp_pe()
Returns processor ID.
Definition: mpp_util.inc:407