| 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 ncdf, not_this => NCDF_SFUN
|
|---|
| 39 |
|
|---|
| 40 | implicit none
|
|---|
| 41 |
|
|---|
| 42 | character(len = *), intent( in) :: name
|
|---|
| 43 | character(len = *), &
|
|---|
| 44 | intent(IN_OR_OUT`') :: values
|
|---|
| 45 | character(len = *), dimension(2), optional :: range
|
|---|
| 46 | character(len = *), optional :: ncfile
|
|---|
| 47 | integer, optional :: ncid
|
|---|
| 48 | integer, optional :: rec
|
|---|
| 49 | integer, dimension(:), optional :: start
|
|---|
| 50 | character(len = *), optional :: units
|
|---|
| 51 |
|
|---|
| 52 | integer :: status, varid, ncid_local, groupid
|
|---|
| 53 | integer :: ndims, dimrec, i
|
|---|
| 54 | integer, dimension(NF90_MAX_VAR_DIMS) :: strt, cnts, dimids
|
|---|
| 55 | character(len = 1), dimension(len(values)+1) &
|
|---|
| 56 | :: rvalues
|
|---|
| 57 | character(len(values)) :: svalues
|
|---|
| 58 | character(len = NF90_MAX_NAME) :: vname
|
|---|
| 59 | logical :: havegroup
|
|---|
| 60 |
|
|---|
| 61 | ifelse(PUTORGET,`put',`
|
|---|
| 62 | ! pgf95 bugfix
|
|---|
| 63 | ! -----------
|
|---|
| 64 |
|
|---|
| 65 | svalues = values
|
|---|
| 66 | '
|
|---|
| 67 | )dnl
|
|---|
| 68 | ifelse(PUTORGET,`get',`
|
|---|
| 69 | ! pgf95 bugfix
|
|---|
| 70 | ! -----------
|
|---|
| 71 |
|
|---|
| 72 | svalues(:) = " "
|
|---|
| 73 | '
|
|---|
| 74 | )dnl
|
|---|
| 75 |
|
|---|
| 76 | ! See if this is the current netcdf
|
|---|
| 77 | ! ---------------------------------
|
|---|
| 78 |
|
|---|
| 79 | if (present(ncfile)) then
|
|---|
| 80 | if (ncfile == ncdf_ncname) then
|
|---|
| 81 | ncid_local = ncdf_ncid
|
|---|
| 82 | else
|
|---|
| 83 | status = nf90_open(ncfile, nf90_share, ncid)
|
|---|
| 84 | if (status /= nf90_noerr) call ncdf_error_handler(status)
|
|---|
| 85 | endif
|
|---|
| 86 |
|
|---|
| 87 | else if (present(ncid)) then
|
|---|
| 88 |
|
|---|
| 89 | ncid_local = ncid
|
|---|
| 90 |
|
|---|
| 91 | else
|
|---|
| 92 |
|
|---|
| 93 | ncid_local = ncdf_ncid
|
|---|
| 94 |
|
|---|
| 95 | endif
|
|---|
| 96 |
|
|---|
| 97 | ! Get group id if necessary and replace ncid_local with it
|
|---|
| 98 | ! --------------------------------------------------------
|
|---|
| 99 |
|
|---|
| 100 | status = ncdf_getgroupid(ncid_local, name, vname, groupid, havegroup)
|
|---|
| 101 | if (status /= nf90_noerr) then
|
|---|
| 102 | WRITE ( *, FMT="(A)" ) "ERROR: Group ID not found: "// name
|
|---|
| 103 | call ncdf_error_handler(status)
|
|---|
| 104 | return
|
|---|
| 105 | endif
|
|---|
| 106 | ncid_local = groupid
|
|---|
| 107 |
|
|---|
| 108 | ! Get variable ID
|
|---|
| 109 | ! ---------------
|
|---|
| 110 |
|
|---|
| 111 | status = nf90_inq_varid(ncid_local, vname, varid)
|
|---|
| 112 | if (status /= nf90_noerr) call ncdf_error_handler(status)
|
|---|
| 113 |
|
|---|
| 114 | ! Obtain some information about the variables dimensionality
|
|---|
| 115 | ! ----------------------------------------------------------
|
|---|
| 116 |
|
|---|
| 117 | status = nf90_inquire_variable(ncid_local, varid, &
|
|---|
| 118 | ndims = ndims, dimids = dimids)
|
|---|
| 119 | if (status /= nf90_noerr) call ncdf_error_handler(status)
|
|---|
| 120 |
|
|---|
| 121 | ! Prepare start and count arrays - these are the defaults
|
|---|
| 122 | ! -------------------------------------------------------
|
|---|
| 123 |
|
|---|
| 124 | strt = 1
|
|---|
| 125 | cnts = 0
|
|---|
| 126 |
|
|---|
| 127 | ! Note: First dim varies along string, second is record...
|
|---|
| 128 | status = nf90_inquire_dimension(ncid_local, dimids(1), len = cnts(1))
|
|---|
| 129 | if (status /= nf90_noerr) call ncdf_error_handler(status)
|
|---|
| 130 | cnts(1) = min(cnts(1), size(rvalues))
|
|---|
| 131 |
|
|---|
| 132 | ! Special cases: record is given...
|
|---|
| 133 | ! ---------------------------------
|
|---|
| 134 |
|
|---|
| 135 | if (present(rec)) then
|
|---|
| 136 |
|
|---|
| 137 | ! ...see if an unlimited (record) dimension is available...
|
|---|
| 138 |
|
|---|
| 139 | status = nf90_inquire(ncid_local, unlimitedDimID = dimrec)
|
|---|
| 140 | if (status /= nf90_noerr .or. dimrec == -1) &
|
|---|
| 141 | call ncdf_error_handler(status)
|
|---|
| 142 |
|
|---|
| 143 | ! ...make sure this is true for the variable in question...
|
|---|
| 144 |
|
|---|
| 145 | i = 1
|
|---|
| 146 | do while (i <= ndims)
|
|---|
| 147 | if (dimids(i) == dimrec) exit
|
|---|
| 148 | i = i + 1
|
|---|
| 149 | enddo
|
|---|
| 150 |
|
|---|
| 151 | ! ...and set the start array.
|
|---|
| 152 |
|
|---|
| 153 | if (i <= ndims) then
|
|---|
| 154 | strt(i) = rec
|
|---|
| 155 | cnts(i) = 1
|
|---|
| 156 | else
|
|---|
| 157 | call ncdf_error_handler(NF90_ENORECVARS)
|
|---|
| 158 | endif
|
|---|
| 159 |
|
|---|
| 160 | ! ...or start...
|
|---|
| 161 | ! --------------
|
|---|
| 162 |
|
|---|
| 163 | else if (present(start)) then
|
|---|
| 164 |
|
|---|
| 165 | strt(1:size(start)) = start
|
|---|
| 166 |
|
|---|
| 167 | endif
|
|---|
| 168 | ifelse(PUTORGET,`put',`
|
|---|
| 169 | ! Copy values
|
|---|
| 170 | ! -----------
|
|---|
| 171 |
|
|---|
| 172 | do i = 1, len_trim(svalues)
|
|---|
| 173 | rvalues(i) = svalues(i:i)
|
|---|
| 174 | enddo
|
|---|
| 175 | if (i <= size(rvalues)) rvalues(i:) = achar(0)
|
|---|
| 176 | '
|
|---|
| 177 | )dnl
|
|---|
| 178 | ifelse(PUTORGET,`get',`
|
|---|
| 179 | ! Initialise values
|
|---|
| 180 | ! -----------------
|
|---|
| 181 |
|
|---|
| 182 | rvalues(:) = achar(0)
|
|---|
| 183 | '
|
|---|
| 184 | )dnl
|
|---|
| 185 | ! Read/write values
|
|---|
| 186 | ! -----------------
|
|---|
| 187 |
|
|---|
| 188 | status = `nf90_'PUTORGET`_var(ncid_local, varid, rvalues, start = strt, count = cnts)'
|
|---|
| 189 | if (status /= nf90_noerr) call ncdf_error_handler(status)
|
|---|
| 190 |
|
|---|
| 191 | ifelse(PUTORGET,`get',`
|
|---|
| 192 | ! Copy values
|
|---|
| 193 | ! -----------
|
|---|
| 194 |
|
|---|
| 195 | do i = 1, len(values)
|
|---|
| 196 | if (iachar(rvalues(i)) > 0) then
|
|---|
| 197 | values(i:i) = rvalues(i)
|
|---|
| 198 | else
|
|---|
| 199 | values(i:i) = " "
|
|---|
| 200 | endif
|
|---|
| 201 | enddo
|
|---|
| 202 |
|
|---|
| 203 | ! Add to counter of number of variables read
|
|---|
| 204 | ! ------------------------------------------
|
|---|
| 205 |
|
|---|
| 206 | if (.not. (havegroup)) ncdf_read(varid) = .true.
|
|---|
| 207 | '
|
|---|
| 208 | )dnl
|
|---|
| 209 |
|
|---|
| 210 | ! Provide formal units
|
|---|
| 211 | ! --------------------
|
|---|
| 212 |
|
|---|
| 213 | if (present(units)) then
|
|---|
| 214 | units = ' '
|
|---|
| 215 | endif
|
|---|
| 216 |
|
|---|
| 217 | ! Dummy lines to avoid warnings - ignore
|
|---|
| 218 | ! --------------------------------------
|
|---|
| 219 | if (present(range)) then
|
|---|
| 220 | continue
|
|---|
| 221 | endif
|
|---|
| 222 |
|
|---|
| 223 | end subroutine NCDF_SFUN
|
|---|