The function NF90_GET_VAR gets one or more data values from a netCDF variable of an open netCDF dataset that is in data mode. Required inputs are the netCDF ID, the variable ID, and a specification for the data values into which the data will be read. Optional inputs may indicate the starting position of the data values in the netCDF variable (argument start), the sampling frequency with which data values are read from the netCDF variable (argument stride), and a mapping between the dimensions of the data array and the netCDF variable (argument map). The values to be read are associated with the netCDF variable by assuming that the first dimension of the netCDF variable varies fastest in the Fortran 90 interface. Data values are converted from the external type of the variable, if necessary.
Take care when using the simplest forms of this interface with record variables when you don't specify how many records are to be read. If you try to read all the values of a record variable into an array but there are more records in the file than you assume, more data will be read than you expect, which may cause a segmentation violation.
function nf90_get_var(ncid, varid, values, start, count, stride, map) integer, intent( in) :: ncid, varid any valid type, scalar or array of any rank, & intent(out) :: values integer, dimension(:), optional, intent( in) :: start, count, stride, map integer :: nf90_get_var
ncid
varid
values
start
By default, start(:) = 1.
count
By default, count(:numDims) = shape(values) and
count(numDims + 1:) = 1, where numDims = size(shape(values)).
stride
By default, stride(:) = 1.
map
By default, edgeLengths = shape(values), and map = (/ 1, (product(edgeLengths(:i)), i = 1, size(edgeLengths) - 1) /), that is, there is no mapping.
Use of Fortran 90 intrinsic functions (including reshape, transpose, and spread) may let you avoid using this argument.
NF90_GET_VAR returns the value NF90_NOERR if no errors occurred. Otherwise, the returned status indicates an error. Possible causes of errors include:
(As noted above, another possible source of error is using this interface to read all the values of a record variable without specifying the number of records. If there are more records in the file than you assume, more data will be read than you expect!)
Here is an example using NF90_GET_VAR to read the (4,3,2) element of the variable named rh from an existing netCDF dataset named foo.nc. For simplicity in this example, we assume that we know that rh is dimensioned with lon, lat, and time, so we want to read the value of rh that corresponds to the fourth lon value, the third lat value, and the second time value:
use netcdf implicit none integer :: ncId, rhVarId, status real :: rhValue ... status = nf90_open("foo.nc", nf90_NoWrite, ncid) if(status /= nf90_NoErr) call handle_err(status) - status = nf90_inq_varid(ncid, "rh", rhVarId) if(status /= nf90_NoErr) call handle_err(status) status = nf90_get_var(ncid, rhVarId, rhValue, start = (/ 4, 3, 2 /) ) if(status /= nf90_NoErr) call handle_err(status)
In this example we use NF90_GET_VAR to read all the values of the variable named rh from an existing netCDF dataset named foo.nc. We assume that we know that rh is dimensioned with lon, lat, and time. In this example we query the netCDF file to discover the lengths of the dimensions, then allocate a Fortran 90 array the same shape as the netCDF variable.
use netcdf implicit none integer :: ncId, rhVarId, & lonDimID, latDimId, timeDimId, & numLons, numLats, numTimes, & status integer, dimension(nf90_max_var_dims) :: dimIDs real, dimension(:, :, :), allocatable :: rhValues ... status = nf90_open("foo.nc", nf90_NoWrite, ncid) if(status /= nf90_NoErr) call handle_err(status) ... status = nf90_inq_varid(ncid, "rh", rhVarId) if(status /= nf90_NoErr) call handle_err(status) ! How big is the netCDF variable, that is, what are the lengths of ! its constituent dimensions? status = nf90_inquire_variable(ncid, rhVarId, dimids = dimIDs) if(status /= nf90_NoErr) call handle_err(status) status = nf90_inquire_dimension(ncid, dimIDs(1), len = numLons) if(status /= nf90_NoErr) call handle_err(status) status = nf90_inquire_dimension(ncid, dimIDs(2), len = numLats) if(status /= nf90_NoErr) call handle_err(status) status = nf90_inquire_dimension(ncid, dimIDs(3), len = numTimes) if(status /= nf90_NoErr) call handle_err(status) allocate(rhValues(numLons, numLats, numTimes)) ... status = nf90_get_var(ncid, rhVarId, rhValues) if(status /= nf90_NoErr) call handle_err(status)
Here is an example using NF90_GET_VAR to read a section of the variable named rh from an existing netCDF dataset named foo.nc. For simplicity in this example, we assume that we know that rh is dimensioned with lon, lat, and time, that there are ten lon values, five lat values, and three time values, and that we want to replace all the values at the last time.
use netcdf implicit none integer :: ncId, rhVarId, status integer, parameter :: numLons = 10, numLats = 5, numTimes = 3 real, dimension(numLons, numLats, numTimes) & :: rhValues ... status = nf90_open("foo.nc", nf90_NoWrite, ncid) if(status /= nf90_NoErr) call handle_err(status) ... status = nf90_inq_varid(ncid, "rh", rhVarId) if(status /= nf90_NoErr) call handle_err(status) !Read the values at the last time by passing an array section status = nf90_get_var(ncid, rhVarId, rhValues(:, :, 3), & start = (/ 1, 1, numTimes /), & count = (/ numLats, numLons, 1 /)) if(status /= nf90_NoErr) call handle_err(status)
Here is an example of using NF_GET_VAR to read every other point of a netCDF variable named rh having dimensions (6, 4).
use netcdf implicit none integer :: ncId, rhVarId, status integer, parameter :: numLons = 6, numLats = 4 real, dimension(numLons, numLats) & :: rhValues ... status = nf90_open("foo.nc", nf90_NoWrite, ncid) if(status /= nf90_NoErr) call handle_err(status) ... status = nf90_inq_varid(ncid, "rh", rhVarId) if(status /= nf90_NoErr) call handle_err(status) ... ! Read every other value into an array section status = nf90_get_var(ncid, rhVarId, rhValues(::2, ::2) & stride = (/ 2, 2 /)) if(status /= nf90_NoErr) call handle_err(status)
The following map vector shows the default mapping between a 2x3x4 netCDF variable and an internal array of the same shape:
real, dimension(2, 3, 4):: a ! same shape as netCDF variable integer, dimension(3) :: map = (/ 1, 2, 6 /) ! netCDF dimension inter-element distance ! ---------------- ---------------------- ! most rapidly varying 1 ! intermediate 2 (= map(1)*2) ! most slowly varying 6 (= map(2)*3)
Using the map vector above obtains the same result as simply not passing a map vector at all.
Here is an example of using nf90_get_var to read a netCDF variable named rh whose dimensions are the transpose of the Fortran 90 array:
use netcdf implicit none integer :: ncId, rhVarId, status integer, parameter :: numLons = 6, numLats = 4 real, dimension(numLons, numLats) :: rhValues ! netCDF variable has dimensions (numLats, numLons) ... status = nf90_open("foo.nc", nf90_NoWrite, ncid) if(status /= nf90_NoErr) call handle_err(status) ... status = nf90_inq_varid(ncid, "rh", rhVarId) if(status /= nf90_NoErr) call handle_err(status) ... ! Read transposed values: map vector would be (/ 1, numLats /) for ! no transposition status = nf90_get_var(ncid, rhVarId,rhValues, map = (/ numLons, 1 /)) if(status /= nf90_NoErr) call handle_err(status)
The same effect can be obtained more simply, though using more memory, using Fortran 90 intrinsic functions:
use netcdf implicit none integer :: ncId, rhVarId, status integer, parameter :: numLons = 6, numLats = 4 real, dimension(numLons, numLats) :: rhValues ! netCDF variable has dimensions (numLats, numLons) real, dimension(numLons, numLats) :: tempValues ... status = nf90_open("foo.nc", nf90_NoWrite, ncid) if(status /= nf90_NoErr) call handle_err(status) ... status = nf90_inq_varid(ncid, "rh", rhVarId) if(status /= nf90_NoErr) call handle_err(status) ... status = nf90_get_var(ncid, rhVarId, tempValues)) if(status /= nf90_NoErr) call handle_err(status) rhValues(:, :) = transpose(tempValues)