| 1 | dnl $Id :$
|
|---|
| 2 | dnl
|
|---|
| 3 | dnl This file is used for the automatic generation of the file
|
|---|
| 4 | dnl ncdf_getvar.f90.
|
|---|
| 5 | dnl
|
|---|
| 6 | dnl AUTHOR
|
|---|
| 7 | dnl C. Marquardt, Darmstadt, Germany <christian@marquardt.sc>
|
|---|
| 8 | dnl
|
|---|
| 9 | dnl COPYRIGHT
|
|---|
| 10 | dnl
|
|---|
| 11 | dnl Copyright (c) 2005 Christian Marquardt <christian@marquardt.sc>
|
|---|
| 12 | dnl
|
|---|
| 13 | dnl All rights reserved.
|
|---|
| 14 | dnl
|
|---|
| 15 | dnl Permission is hereby granted, free of charge, to any person obtaining
|
|---|
| 16 | dnl a copy of this software and associated documentation files (the
|
|---|
| 17 | dnl "Software"), to deal in the Software without restriction, including
|
|---|
| 18 | dnl without limitation the rights to use, copy, modify, merge, publish,
|
|---|
| 19 | dnl distribute, sublicense, and/or sell copies of the Software, and to
|
|---|
| 20 | dnl permit persons to whom the Software is furnished to do so, subject to
|
|---|
| 21 | dnl the following conditions:
|
|---|
| 22 | dnl
|
|---|
| 23 | dnl The above copyright notice and this permission notice shall be
|
|---|
| 24 | dnl included in all copies or substantial portions of the Software as well
|
|---|
| 25 | dnl as in supporting documentation.
|
|---|
| 26 | dnl
|
|---|
| 27 | dnl THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND,
|
|---|
| 28 | dnl EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF
|
|---|
| 29 | dnl MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND
|
|---|
| 30 | dnl NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE
|
|---|
| 31 | dnl LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION
|
|---|
| 32 | dnl OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION
|
|---|
| 33 | dnl WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
|
|---|
| 34 | dnl
|
|---|
| 35 | subroutine NCDF_SFUN (name, values, ncfile, ncid, rec, start, units, range)
|
|---|
| 36 |
|
|---|
| 37 | use typeSizes
|
|---|
| 38 | use unitconvert, only: ut_convert
|
|---|
| 39 | use ncdf, not_this => NCDF_SFUN
|
|---|
| 40 |
|
|---|
| 41 | implicit none
|
|---|
| 42 |
|
|---|
| 43 | character(len = *), intent( in) :: name
|
|---|
| 44 | TYPE, &
|
|---|
| 45 | intent(IN_OR_OUT`') :: values
|
|---|
| 46 | TYPE, dimension(2), optional :: range
|
|---|
| 47 | character(len = *), optional :: ncfile
|
|---|
| 48 | integer, optional :: ncid
|
|---|
| 49 | integer, optional :: rec
|
|---|
| 50 | integer, dimension(:), optional :: start
|
|---|
| 51 | character(len = *), optional :: units
|
|---|
| 52 |
|
|---|
| 53 | integer :: status, varid, groupid, ncid_local
|
|---|
| 54 | integer :: ndims, dimrec, i
|
|---|
| 55 | integer, dimension(NF90_MAX_VAR_DIMS) :: strt, dimids
|
|---|
| 56 | character(len = NF90_MAX_NAME) :: ncdf_units
|
|---|
| 57 | TYPE &
|
|---|
| 58 | :: rvalues
|
|---|
| 59 | TYPE, dimension(2) :: rrange
|
|---|
| 60 | TYPE &
|
|---|
| 61 | :: scale_factor, add_offset
|
|---|
| 62 | logical :: have_scale, have_offset
|
|---|
| 63 | logical :: have_range
|
|---|
| 64 | character(len = NF90_MAX_NAME) :: vname
|
|---|
| 65 | logical :: havegroup
|
|---|
| 66 |
|
|---|
| 67 | ! See if this is the current netcdf
|
|---|
| 68 | ! ---------------------------------
|
|---|
| 69 |
|
|---|
| 70 | if (present(ncfile)) then
|
|---|
| 71 | if (ncfile == ncdf_ncname) then
|
|---|
| 72 | ncid_local = ncdf_ncid
|
|---|
| 73 | else
|
|---|
| 74 | status = nf90_open(ncfile, nf90_share, ncid)
|
|---|
| 75 | if (status /= nf90_noerr) call ncdf_error_handler(status)
|
|---|
| 76 | endif
|
|---|
| 77 |
|
|---|
| 78 | else if (present(ncid)) then
|
|---|
| 79 |
|
|---|
| 80 | ncid_local = ncid
|
|---|
| 81 |
|
|---|
| 82 | else
|
|---|
| 83 |
|
|---|
| 84 | ncid_local = ncdf_ncid
|
|---|
| 85 |
|
|---|
| 86 | endif
|
|---|
| 87 |
|
|---|
| 88 | ! Get group id if necessary and replace ncid_local with it
|
|---|
| 89 | ! --------------------------------------------------------
|
|---|
| 90 |
|
|---|
| 91 | status = ncdf_getgroupid(ncid_local, name, vname, groupid, havegroup)
|
|---|
| 92 | if (status /= nf90_noerr) then
|
|---|
| 93 | WRITE ( *, FMT="(A)" ) "ERROR: Group ID not found: "// name
|
|---|
| 94 | call ncdf_error_handler(status)
|
|---|
| 95 | return
|
|---|
| 96 | endif
|
|---|
| 97 | ncid_local = groupid
|
|---|
| 98 |
|
|---|
| 99 | ! Get variable ID
|
|---|
| 100 | ! ---------------
|
|---|
| 101 |
|
|---|
| 102 | status = nf90_inq_varid(ncid_local, vname, varid)
|
|---|
| 103 | if (status /= nf90_noerr) then
|
|---|
| 104 | WRITE ( *, FMT="(A)" ) "ERROR: Variable not found: "// name
|
|---|
| 105 | call ncdf_error_handler(status)
|
|---|
| 106 | return
|
|---|
| 107 | endif
|
|---|
| 108 |
|
|---|
| 109 | ! Obtain some information about the variables dimensionality
|
|---|
| 110 | ! ----------------------------------------------------------
|
|---|
| 111 |
|
|---|
| 112 | status = nf90_inquire_variable(ncid_local, varid, &
|
|---|
| 113 | ndims = ndims, dimids = dimids)
|
|---|
| 114 | if (status /= nf90_noerr) call ncdf_error_handler(status)
|
|---|
| 115 |
|
|---|
| 116 | ! Obtain scaling factors for the variable
|
|---|
| 117 | ! ---------------------------------------
|
|---|
| 118 |
|
|---|
| 119 | status = nf90_get_att(ncid_local, varid, 'scale_factor', scale_factor)
|
|---|
| 120 | if (status == nf90_enotatt) then
|
|---|
| 121 | have_scale = .false.
|
|---|
| 122 | status = nf90_noerr
|
|---|
| 123 | else
|
|---|
| 124 | have_scale = .true.
|
|---|
| 125 | endif
|
|---|
| 126 | if (status /= nf90_noerr) then
|
|---|
| 127 | WRITE ( *, FMT="(A)" ) "ERROR: Attribute scale_factor not found for variable: "// name
|
|---|
| 128 | call ncdf_error_handler(status)
|
|---|
| 129 | endif
|
|---|
| 130 |
|
|---|
| 131 | status = nf90_get_att(ncid_local, varid, 'add_offset', add_offset)
|
|---|
| 132 | if (status == nf90_enotatt) then
|
|---|
| 133 | have_offset = .false.
|
|---|
| 134 | status = nf90_noerr
|
|---|
| 135 | else
|
|---|
| 136 | have_offset = .true.
|
|---|
| 137 | endif
|
|---|
| 138 | if (status /= nf90_noerr) then
|
|---|
| 139 | WRITE ( *, FMT="(A)" ) "ERROR: Attribute add_offset not found for variable: "// name
|
|---|
| 140 | call ncdf_error_handler(status)
|
|---|
| 141 | endif
|
|---|
| 142 | ! Prepare start and count arrays - these are the defaults
|
|---|
| 143 | ! -------------------------------------------------------
|
|---|
| 144 |
|
|---|
| 145 | strt = 1
|
|---|
| 146 |
|
|---|
| 147 | ! Special cases: record is given...
|
|---|
| 148 | ! ---------------------------------
|
|---|
| 149 |
|
|---|
| 150 | if (present(rec)) then
|
|---|
| 151 |
|
|---|
| 152 | ! ...see if an unlimited (record) dimension is available...
|
|---|
| 153 |
|
|---|
| 154 | status = nf90_inquire(ncid_local, unlimitedDimID = dimrec)
|
|---|
| 155 | if (status /= nf90_noerr .or. dimrec == -1) &
|
|---|
| 156 | call ncdf_error_handler(status)
|
|---|
| 157 |
|
|---|
| 158 | ! ...make sure this is true for the variable in question...
|
|---|
| 159 |
|
|---|
| 160 | i = 1
|
|---|
| 161 | do while (i <= ndims)
|
|---|
| 162 | if (dimids(i) == dimrec) exit
|
|---|
| 163 | i = i + 1
|
|---|
| 164 | enddo
|
|---|
| 165 |
|
|---|
| 166 | ! ...and set the start array.
|
|---|
| 167 |
|
|---|
| 168 | if (i <= ndims) then
|
|---|
| 169 | strt(i) = rec
|
|---|
| 170 | else
|
|---|
| 171 | call ncdf_error_handler(NF90_ENORECVARS)
|
|---|
| 172 | endif
|
|---|
| 173 |
|
|---|
| 174 | ! ...or start...
|
|---|
| 175 | ! --------------
|
|---|
| 176 |
|
|---|
| 177 | else if (present(start)) then
|
|---|
| 178 |
|
|---|
| 179 | strt(1:size(start)) = start
|
|---|
| 180 |
|
|---|
| 181 | endif
|
|---|
| 182 |
|
|---|
| 183 | ifelse(PUTORGET,`put',dnl
|
|---|
| 184 | ! Convert units
|
|---|
| 185 | ! -------------
|
|---|
| 186 |
|
|---|
| 187 | if (present(units)) then
|
|---|
| 188 | ncdf_units(:) = ' '
|
|---|
| 189 | status = nf90_get_att(ncid_local, varid, 'units', ncdf_units)
|
|---|
| 190 | if (status /= nf90_noerr) then
|
|---|
| 191 | WRITE ( *, FMT="(A)" ) "ERROR: Attribute units not found for variable: "// name
|
|---|
| 192 | call ncdf_error_handler(status)
|
|---|
| 193 | endif
|
|---|
| 194 | call ut_convert(values, units, rvalues, ncdf_units)
|
|---|
| 195 |
|
|---|
| 196 | ! Only convert valid data
|
|---|
| 197 | status = nf90_get_att(ncid_local, varid, 'valid_range', rrange)
|
|---|
| 198 | if (status /= nf90_noerr) then
|
|---|
| 199 | have_range=.false.
|
|---|
| 200 | else
|
|---|
| 201 | have_range=.true.
|
|---|
| 202 | if(rvalues < rrange(1) .or. rvalues > rrange(2)) rvalues = values
|
|---|
| 203 | endif
|
|---|
| 204 |
|
|---|
| 205 | else
|
|---|
| 206 | rvalues = values
|
|---|
| 207 | endif
|
|---|
| 208 |
|
|---|
| 209 | ! Scale variables if necessary
|
|---|
| 210 | ! ----------------------------
|
|---|
| 211 |
|
|---|
| 212 | if (have_scale .and. have_offset) then
|
|---|
| 213 | rvalues = (rvalues - add_offset) / scale_factor
|
|---|
| 214 | else if (have_scale .and. (.not. have_offset)) then
|
|---|
| 215 | rvalues = rvalues / scale_factor
|
|---|
| 216 | else if ((.not. have_scale) .and. have_offset) then
|
|---|
| 217 | rvalues = rvalues - add_offset
|
|---|
| 218 | endif
|
|---|
| 219 | ,dnl
|
|---|
| 220 | ! Copy data values
|
|---|
| 221 | ! ----------------
|
|---|
| 222 |
|
|---|
| 223 | rvalues = values
|
|---|
| 224 | )
|
|---|
| 225 | ! Read/write values
|
|---|
| 226 | ! -----------------
|
|---|
| 227 |
|
|---|
| 228 | status = `nf90_'PUTORGET`_var(ncid_local, varid, rvalues, start = strt)'
|
|---|
| 229 | if (status /= nf90_noerr) call ncdf_error_handler(status)
|
|---|
| 230 |
|
|---|
| 231 | ifelse(PUTORGET,`get',dnl
|
|---|
| 232 | ! Scale variables if necessary
|
|---|
| 233 | ! ----------------------------
|
|---|
| 234 |
|
|---|
| 235 | if (have_scale .and. have_offset) then
|
|---|
| 236 | rvalues = rvalues * scale_factor + add_offset
|
|---|
| 237 | else if (have_scale .and. (.not. have_offset)) then
|
|---|
| 238 | rvalues = rvalues * scale_factor
|
|---|
| 239 | else if ((.not. have_scale) .and. have_offset) then
|
|---|
| 240 | rvalues = rvalues + add_offset
|
|---|
| 241 | endif
|
|---|
| 242 |
|
|---|
| 243 | ! Obtain range values - use defaults if not present
|
|---|
| 244 | ! -------------------
|
|---|
| 245 |
|
|---|
| 246 | if (present(range))then
|
|---|
| 247 | status = nf90_get_att(ncid_local, varid, 'valid_range', rrange)
|
|---|
| 248 | if (status /= nf90_noerr) then
|
|---|
| 249 | rrange = range
|
|---|
| 250 | have_range = .false.
|
|---|
| 251 | else
|
|---|
| 252 | have_range = .true.
|
|---|
| 253 | endif
|
|---|
| 254 | endif
|
|---|
| 255 |
|
|---|
| 256 | ! Convert units - if present & not the same as target
|
|---|
| 257 | ! -------------
|
|---|
| 258 |
|
|---|
| 259 | if (present(units)) then
|
|---|
| 260 | ncdf_units(:) = ' '
|
|---|
| 261 | status = nf90_get_att(ncid_local, varid, 'units', ncdf_units)
|
|---|
| 262 | if (status /= nf90_noerr) then
|
|---|
| 263 | WRITE ( *, FMT="(A)" ) "ERROR: Attribute units not found for variable: "// name
|
|---|
| 264 | call ncdf_error_handler(status)
|
|---|
| 265 | endif
|
|---|
| 266 | call ut_convert(rvalues, ncdf_units, values, units)
|
|---|
| 267 |
|
|---|
| 268 | if(present(range) .and. have_range) then
|
|---|
| 269 | call ut_convert(rrange, ncdf_units, range, units)
|
|---|
| 270 |
|
|---|
| 271 | ! Only convert valid data
|
|---|
| 272 | if(rvalues < rrange(1) .or. rvalues > rrange(2)) values = rvalues
|
|---|
| 273 | endif
|
|---|
| 274 |
|
|---|
| 275 | else
|
|---|
| 276 | values = rvalues
|
|---|
| 277 | if(present(range)) range = rrange
|
|---|
| 278 | endif
|
|---|
| 279 |
|
|---|
| 280 | ! Add to counter of number of variables read (deactivated for n4)
|
|---|
| 281 | ! ------------------------------------------
|
|---|
| 282 |
|
|---|
| 283 | if (.not. (havegroup)) ncdf_read(varid) = .true.
|
|---|
| 284 |
|
|---|
| 285 | )dnl
|
|---|
| 286 |
|
|---|
| 287 | ! Dummy lines to avoid warnings - ignore
|
|---|
| 288 | ! --------------------------------------
|
|---|
| 289 | if (present(range)) then
|
|---|
| 290 | continue
|
|---|
| 291 | endif
|
|---|
| 292 |
|
|---|
| 293 | end subroutine NCDF_SFUN
|
|---|