Next: NF90_RENAME_ATT, Previous: NF90_GET_ATT, Up: Attributes
The function NF90_COPY_ATT copies an attribute from one open netCDF dataset to another. It can also be used to copy an attribute from one variable to another within the same netCDF dataset.
function nf90_copy_att(ncid_in, varid_in, name, ncid_out, varid_out)
integer, intent( in) :: ncid_in, varid_in
character (len = *), intent( in) :: name
integer, intent( in) :: ncid_out, varid_out
integer :: nf90_copy_att
ncid_invarid_innamencid_outvarid_outNF90_COPY_ATT returns the value NF90_NOERR if no errors occurred. Otherwise, the returned status indicates an error. Possible causes of errors include:
Here is an example using NF90_COPY_ATT to copy the variable attribute units from the variable rh in an existing netCDF dataset named foo.nc to the variable avgrh in another existing netCDF dataset named bar.nc, assuming that the variable avgrh already exists, but does not yet have a units attribute:
use netcdf
implicit none
integer :: ncid1, ncid2, status
integer :: RHVarID, avgRHVarID ! Variable ID
...
status = nf90_open("foo.nc", nf90_nowrite, ncid1)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_open("bar.nc", nf90_write, ncid2)
if (status /= nf90_noerr) call handle_err(status)
...
! Find the IDs of the variables
status = nf90_inq_varid(ncid1, "rh", RHVarID)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_inq_varid(ncid1, "avgrh", avgRHVarID)
if (status /= nf90_noerr) call handle_err(status)
...
status = nf90_redef(ncid2) ! Enter define mode
if (status /= nf90_noerr) call handle_err(status)
! Copy variable attribute from "rh" in file 1 to "avgrh" in file 1
status = nf90_copy_att(ncid1, RHVarID, "units", ncid2, avgRHVarID)
if (status /= nf90_noerr) call handle_err(status)
status = nf90_enddef(ncid2)
if (status /= nf90_noerr) call handle_err(status)