Skip to content
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
3 changes: 2 additions & 1 deletion configure.ac
Original file line number Diff line number Diff line change
Expand Up @@ -683,7 +683,8 @@ AC_SUBST(NF_HDF5_PLUGIN_PATH, [$nc_hdf5_plugin_path])
# to their source directory, or the test program makes an assumption
# about where files live. AC_CONFIG_LINKS provides a mechanism to
# link/copy files if an out-of-source build is happening.
AC_CONFIG_LINKS([nf_test/fills.nc:nf_test/ref_fills.nc])
AC_CONFIG_LINKS([nf_test/fills.nc:nf_test/ref_fills.nc
nf03_test4/att.nc:nf03_test4/ref_att.nc])

if test "x$nc_has_parallel" = xyes; then
AC_CONFIG_FILES([nf03_test4/run_f90_par_test.sh], [chmod ugo+x nf03_test4/run_f90_par_test.sh])
Expand Down
48 changes: 44 additions & 4 deletions fortran/netcdf_attributes.F90
Original file line number Diff line number Diff line change
Expand Up @@ -67,13 +67,53 @@ function nf90_put_att_text(ncid, varid, name, values)
end function nf90_put_att_text
! -------
function nf90_get_att_text(ncid, varid, name, values)
use, intrinsic :: iso_c_binding, only: c_ptr, c_size_t, c_f_pointer, c_int
implicit none
integer, intent( in) :: ncid, varid
character(len = *), intent( in) :: name
character(len = *), intent(out) :: values
integer :: nf90_get_att_text

values = ' ' !! make sure result will be blank padded
nf90_get_att_text = nf_get_att_text(ncid, varid, name, values)
interface
integer(c_int) function nc_get_att_string(ncid, varid, name, pp) bind(c)
use iso_c_binding , only : c_int , c_char , c_ptr
integer(c_int) , value :: ncid , varid
character(kind=c_char) , intent(in) :: name
type(c_ptr), intent(out) :: pp
end function nc_get_att_string
end interface
interface
integer(c_size_t) function strlen(cs) bind(c, name='strlen')
use, intrinsic :: iso_c_binding , only : c_size_t , c_ptr
implicit none
type(c_ptr), intent(in), value :: cs
end function strlen
end interface
integer :: xtype , nlen , attid , i
integer(c_int) :: c_ncid , c_varid , c_status , c_nlen
type(c_ptr) :: c_str
character(len_trim(name)+1) :: c_aname
character , pointer :: f_str(:)
nf90_get_att_text = nf90_inquire_attribute(ncid, varid, name, &
xtype, nlen, attid)
if ( nf90_get_att_text == nf90_noerr ) then
if ( xtype == nf90_string .and. nlen == 1 ) then
c_ncid = ncid
c_varid = varid - 1
c_aname = name//char(0)
c_status = nc_get_att_string(c_ncid, c_varid, c_aname, c_str)
nf90_get_att_text = c_status
if ( nf90_get_att_text == nf90_noerr ) then
call c_f_pointer(c_str,f_str,[strlen(c_str)])
values = adjustl("")
do i = 1, size(f_str)
values(i:i) = f_str(i)
end do
end if
else
values = ' ' !! make sure result will be blank padded
nf90_get_att_text = nf_get_att_text(ncid,varid,name,values)
end if
end if
end function nf90_get_att_text
! -------
! Integer attributes
Expand Down Expand Up @@ -183,7 +223,7 @@ function nf90_get_att_FourByteInt(ncid, varid, name, values)
character(len = *), intent( in) :: name
integer (kind = FourByteInt), dimension(:), intent(out) :: values
integer :: nf90_get_att_FourByteInt

integer, dimension(size(values)) :: defaultInteger

nf90_get_att_FourByteInt = nf_get_att_int(ncid, varid, name, defaultInteger)
Expand Down
10 changes: 9 additions & 1 deletion nf03_test4/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -7,7 +7,8 @@ if (USE_NETCDF4)
# Add these netCDF-4 test programs.
SET(nc4_check_PROGRAMS f90tst_vars f90tst_vars_vlen f90tst_grps
f90tst_fill f90tst_fill2 f90tst_vars3 f90tst_vars4 f90tst_vars2
f90tst_path f90tst_rengrps f90tst_nc4 f90tst_types f90tst_types2 f90tst_zstandard)
f90tst_path f90tst_rengrps f90tst_nc4 f90tst_types f90tst_types2 f90tst_zstandard
f90tst_att)
# Only run this test if quantize is present (i.e. netcdf-c is 4.9.0 or greater).
IF (HAVE_QUANTIZE)
SET(nc4_check_PROGRAMS ${nc4_check_PROGRAMS} f90tst_vars5)
Expand All @@ -30,6 +31,13 @@ if (USE_NETCDF4)
SET(f90tst_types_SOURCES f90tst_types.F90)
SET(f90tst_types2_SOURCES f90tst_types2.F90)
SET(f90tst_zstandard_SOURCES f90tst_zstandard.F90)
SET(f90tst_att_SOURCES f90tst_att.F90)

# Need a copy of ref_atts.nc for f90tst_att
file(COPY ${CMAKE_CURRENT_SOURCE_DIR}/ref_att.nc
DESTINATION ${CMAKE_CURRENT_BINARY_DIR})
file(RENAME ${CMAKE_CURRENT_BINARY_DIR}/ref_att.nc
${CMAKE_CURRENT_BINARY_DIR}/att.nc)

# Test parallel I/O.
IF (TEST_PARALLEL)
Expand Down
12 changes: 9 additions & 3 deletions nf03_test4/Makefile.am
Original file line number Diff line number Diff line change
Expand Up @@ -15,7 +15,8 @@ LDADD = ${top_builddir}/fortran/libnetcdff.la
# tst_f90_nc4
NC4_F90_TESTS = f90tst_vars f90tst_vars_vlen f90tst_grps f90tst_fill \
f90tst_fill2 f90tst_vars2 f90tst_vars3 f90tst_vars4 f90tst_path \
f90tst_rengrps f90tst_nc4 f90tst_types f90tst_types2 f90tst_zstandard
f90tst_rengrps f90tst_nc4 f90tst_types f90tst_types2 f90tst_zstandard \
f90tst_att

if TEST_QUANTIZE
NC4_F90_TESTS += f90tst_vars5
Expand All @@ -33,6 +34,7 @@ f90tst_vars2_SOURCES = f90tst_vars2.F90
f90tst_vars3_SOURCES = f90tst_vars3.F90
f90tst_vars4_SOURCES = f90tst_vars4.F90
f90tst_vars5_SOURCES = f90tst_vars5.F90
f90tst_att_SOURCES = f90tst_att.F90
f90tst_path_SOURCES = f90tst_path.F90
f90tst_rengrps_SOURCES = f90tst_rengrps.F90
f90tst_nc4_SOURCES = f90tst_nc4.F90
Expand Down Expand Up @@ -81,7 +83,11 @@ endif # TEST_SZIP_WRITE
@VALGRIND_CHECK_RULES@

EXTRA_DIST = CMakeLists.txt f90tst_flarge.F90 \
run_f90_par_test.sh.in
run_f90_par_test.sh.in \
ref_att.cdl ref_att.nc

# Cleaning up files created during the testing.
CLEANFILES = f90tst_*.nc fort.*
CLEANFILES = f90tst_*.nc fort.*

# This file gets copied to build dir by configure.
DISTCLEANFILES = att.nc
75 changes: 75 additions & 0 deletions nf03_test4/f90tst_att.F90
Original file line number Diff line number Diff line change
@@ -0,0 +1,75 @@
! This program tests reading NF90_STRING & NF90_CHAR attributes from fortran.
! The ref input data structure is shown in ref_att.cdl
! The ref input data can be generated by running "ncgen -o ref_att.nc ref_att.cdl"
! For the input of this program, rename ref input data as "att.nc"
!
! Cheng Da 2025/08/30
!
program f90tst_att
use netcdf, only: NF90_NOWRITE,NF90_GLOBAL,nf90_open,nf90_close,&
nf90_inq_ncid, nf90_inq_varid,nf90_get_att
implicit none

integer :: fid, gid, vid
integer :: ierr
character(len=256) :: fname, str

fname = "att.nc"
call check(NF90_Open(TRIM(fname), NF90_NOWRITE, fid))
! read global att
str=""; call check(nf90_get_att(fid, NF90_GLOBAL, "agency", str)) !NF90_CHAR att
call verify_results(trim(str), "fake_agency", "get wrong global attr: agency", 10)
str=""; call check(nf90_get_att(fid, NF90_GLOBAL, "author", str)) !NF90_STRING att
call verify_results(trim(str), "fake_author", "get wrong global attr: author", 11)

! read var att
call check(nf90_inq_varid(fid, "obstype", vid))
str=""; call check(nf90_get_att(fid, vid, "long_name", str)) !NF90_STRING att
call verify_results(trim(str), "observation type", "get wrong obstype attr: long_name", 12)
str=""; call check(nf90_get_att(fid, vid, "units", str)) !NF90_CHAR att
call verify_results(trim(str), "N/A", "get wrong obstype attr: units", 13)

! read var atts under group MetaData
call check(nf90_inq_ncid(fid, "MetaData", gid))
call check(nf90_inq_varid(gid,"latitude", vid))
str=""; call check(nf90_get_att(gid, vid, "units", str)) ! NF90_STRING att
call verify_results(trim(str), "degrees_north", "get wrong MetaData/latitude attr: units", 14)

! read var atts group yobs
call check(nf90_inq_ncid(fid, "yobs", gid))
call check(nf90_inq_varid(gid,"sst", vid))
str=""; call check(nf90_get_att(gid, vid, "units", str)) ! NF90_STRING att
call verify_results(trim(str), "C", "get wrong yobs/sst attr: units", 15)

! read var atts group Hx
call check(nf90_inq_ncid(fid, "Hx", gid))
call check(nf90_inq_varid(gid,"sst", vid))
str=""; call check(nf90_get_att(gid, vid, "units", str)) ! NF90_CHAR att
call verify_results(trim(str), "C", "get wrong Hx/sst attr: units", 16)
call check(nf90_close(fid))
print *,"*** SUCCESS!"

contains
subroutine check(errcode)
use netcdf, only: NF90_NOERR, nf90_strerror
implicit none
integer, intent(in) :: errcode
if(errcode /= NF90_NOERR) then
print*, "Error: ", trim(nf90_strerror(errcode))
stop (1)
endif
end subroutine check

subroutine verify_results(input, ref, errmsg, errcode)
implicit none
character(*),intent(in) :: input, ref, errmsg
integer,intent(in) :: errcode
if (trim(input) .ne. trim(ref)) then
print*, "Error: "//trim(errmsg)
print*, " expect to get:"//trim(ref)
print*, " actually get:"//trim(input)
stop (errcode)
end if
end subroutine verify_results

end program f90tst_att
32 changes: 32 additions & 0 deletions nf03_test4/ref_att.cdl
Original file line number Diff line number Diff line change
@@ -0,0 +1,32 @@
netcdf ref_att {

dimensions:
Location = UNLIMITED ; // (0 currently)

// global attributes:
string :author = "fake_author" ;
:agency = "fake_agency" ;

variables:
int obstype(Location) ;
string obstype:long_name = "observation type" ;
obstype:units = "N/A";

group: MetaData {
variables:
float latitude(Location) ;
string latitude:units = "degrees_north" ;
}

group: yobs {
variables:
float sst(Location) ;
string sst:units = "C" ;
}

group: Hx {
variables:
float sst(Location) ;
sst:units = "C" ;
}
}
Binary file added nf03_test4/ref_att.nc
Binary file not shown.