Skip to content

Commit

Permalink
Fix for writing single values into 1D array
Browse files Browse the repository at this point in the history
For example, appending a value to an unlimited variable
  • Loading branch information
ZedThree committed Feb 28, 2022
1 parent 911ebbf commit 13f188d
Show file tree
Hide file tree
Showing 3 changed files with 100 additions and 10 deletions.
37 changes: 32 additions & 5 deletions src/neasyf.f90
Original file line number Diff line number Diff line change
Expand Up @@ -868,16 +868,20 @@ subroutine neasyf_dim(parent_id, name, values, dim_size, dimid, varid, units, lo
end if
end subroutine neasyf_dim

subroutine neasyf_write_scalar(parent_id, name, values, units, long_name, start)
subroutine neasyf_write_scalar(parent_id, name, values, dim_ids, dim_names, units, long_name, start)
use, intrinsic :: iso_fortran_env, only : int8, int16, int32, real32, real64
use netcdf, only : nf90_inq_varid, nf90_def_var, nf90_put_var, nf90_put_att, &
NF90_NOERR, NF90_ENOTVAR, nf90_def_dim, NF90_CHAR
NF90_NOERR, NF90_ENOTVAR, NF90_EDIMMETA, nf90_def_dim, NF90_CHAR, nf90_inq_dimid
!> Name of the variable
character(len=*), intent(in) :: name
!> NetCDF ID of the parent group/file
integer, intent(in) :: parent_id
!> Value of the integer to write
class(*), intent(in) :: values
!> Array of dimension IDs
integer, dimension(1), optional, intent(in) :: dim_ids
!> Array of dimension names
character(len=*), dimension(1), optional, intent(in) :: dim_names
!> Units of coordinate
character(len=*), optional, intent(in) :: units
!> Long descriptive name
Expand All @@ -888,6 +892,7 @@ subroutine neasyf_write_scalar(parent_id, name, values, units, long_name, start)
integer(nf_kind) :: nf_type
integer :: status
integer :: var_id, dim_id
integer, dimension(1) :: local_dim_ids

status = nf90_inq_varid(parent_id, name, var_id)
! Variable doesn't exist, so let's create it
Expand All @@ -899,11 +904,33 @@ subroutine neasyf_write_scalar(parent_id, name, values, units, long_name, start)

status = nf90_def_var(parent_id, name, NF90_CHAR, [dim_id], var_id)
class default
nf_type = neasyf_type(values)
! TODO: check if nf_type indicates a derived type
status = nf90_def_var(parent_id, name, nf_type, var_id)
call neasyf_error(status, var=name, varid=var_id, message="defining variable")
nf_type = neasyf_type(values)

if (present(start)) then
! If we've been passed `start`, then we must be writing into
! a 1D array, in which case we need to know the associated
! dimension in order to create the variable
if (.not. (present(dim_ids) .or. present(dim_names))) then
call neasyf_error(NF90_EDIMMETA, var=name, &
message="neasyf_write: one of 'dim_ids' or 'dim_names' must be present if 1D variable doesn't already exist")
end if

! Either use the passed in array of dimension IDs, or look them up from the array of names
if (present(dim_ids)) then
local_dim_ids = dim_ids
else
status = nf90_inq_dimid(parent_id, trim(dim_names(1)), local_dim_ids(1))
call neasyf_error(status, var=name, dim=trim(dim_names(1)))
end if
! Create a 1D variable, although we're only going to write a single value
status = nf90_def_var(parent_id, name, nf_type, local_dim_ids, var_id)
else
! Create a scalar variable
status = nf90_def_var(parent_id, name, nf_type, var_id)
end if
end select
call neasyf_error(status, var=name, varid=var_id, message="defining variable")

if (present(units)) then
status = nf90_put_att(parent_id, var_id, "units", units)
Expand Down
37 changes: 32 additions & 5 deletions src/neasyf.in.f90
Original file line number Diff line number Diff line change
Expand Up @@ -428,16 +428,20 @@ subroutine neasyf_dim(parent_id, name, values, dim_size, dimid, varid, units, lo
end if
end subroutine neasyf_dim

subroutine neasyf_write_scalar(parent_id, name, values, units, long_name, start)
subroutine neasyf_write_scalar(parent_id, name, values, dim_ids, dim_names, units, long_name, start)
use, intrinsic :: iso_fortran_env, only : int8, int16, int32, real32, real64
use netcdf, only : nf90_inq_varid, nf90_def_var, nf90_put_var, nf90_put_att, &
NF90_NOERR, NF90_ENOTVAR, nf90_def_dim, NF90_CHAR
NF90_NOERR, NF90_ENOTVAR, NF90_EDIMMETA, nf90_def_dim, NF90_CHAR, nf90_inq_dimid
!> Name of the variable
character(len=*), intent(in) :: name
!> NetCDF ID of the parent group/file
integer, intent(in) :: parent_id
!> Value of the integer to write
class(*), intent(in) :: values
!> Array of dimension IDs
integer, dimension(1), optional, intent(in) :: dim_ids
!> Array of dimension names
character(len=*), dimension(1), optional, intent(in) :: dim_names
!> Units of coordinate
character(len=*), optional, intent(in) :: units
!> Long descriptive name
Expand All @@ -448,6 +452,7 @@ subroutine neasyf_write_scalar(parent_id, name, values, units, long_name, start)
integer(nf_kind) :: nf_type
integer :: status
integer :: var_id, dim_id
integer, dimension(1) :: local_dim_ids

status = nf90_inq_varid(parent_id, name, var_id)
! Variable doesn't exist, so let's create it
Expand All @@ -459,11 +464,33 @@ subroutine neasyf_write_scalar(parent_id, name, values, units, long_name, start)

status = nf90_def_var(parent_id, name, NF90_CHAR, [dim_id], var_id)
class default
nf_type = neasyf_type(values)
! TODO: check if nf_type indicates a derived type
status = nf90_def_var(parent_id, name, nf_type, var_id)
call neasyf_error(status, var=name, varid=var_id, message="defining variable")
nf_type = neasyf_type(values)

if (present(start)) then
! If we've been passed `start`, then we must be writing into
! a 1D array, in which case we need to know the associated
! dimension in order to create the variable
if (.not. (present(dim_ids) .or. present(dim_names))) then
call neasyf_error(NF90_EDIMMETA, var=name, &
message="neasyf_write: one of 'dim_ids' or 'dim_names' must be present if 1D variable doesn't already exist")
end if

! Either use the passed in array of dimension IDs, or look them up from the array of names
if (present(dim_ids)) then
local_dim_ids = dim_ids
else
status = nf90_inq_dimid(parent_id, trim(dim_names(1)), local_dim_ids(1))
call neasyf_error(status, var=name, dim=trim(dim_names(1)))
end if
! Create a 1D variable, although we're only going to write a single value
status = nf90_def_var(parent_id, name, nf_type, local_dim_ids, var_id)
else
! Create a scalar variable
status = nf90_def_var(parent_id, name, nf_type, var_id)
end if
end select
call neasyf_error(status, var=name, varid=var_id, message="defining variable")

if (present(units)) then
status = nf90_put_att(parent_id, var_id, "units", units)
Expand Down
36 changes: 36 additions & 0 deletions tests/test_write.pf
Original file line number Diff line number Diff line change
Expand Up @@ -30,6 +30,42 @@ contains
close(unit, status="delete")
end subroutine teardown

@test
subroutine test_write_scalar
integer, parameter :: read_data = 5
integer :: expected_data
integer :: var_id

call neasyf_write(file_id, "data", read_data)
call neasyf_close(file_id)

call neasyf_error(nf90_open(temp_file, NF90_NOWRITE, file_id))
call neasyf_error(nf90_inq_varid(file_id, "data", var_id))
call neasyf_error(nf90_get_var(file_id, var_id, expected_data))
call neasyf_close(file_id)

@assertEqual(expected_data, read_data)
end subroutine test_write_scalar

@test
subroutine test_write_scalar_unlimited
integer, dimension(2), parameter :: read_data = [5, 6]
integer, dimension(2) :: expected_data
integer :: var_id

call neasyf_dim(file_id, "t", unlimited=.true.)
call neasyf_write(file_id, "t", read_data(1), dim_names=["t"], start=[1])
call neasyf_write(file_id, "t", read_data(2), dim_names=["t"], start=[2])
call neasyf_close(file_id)

call neasyf_error(nf90_open(temp_file, NF90_NOWRITE, file_id))
call neasyf_error(nf90_inq_varid(file_id, "t", var_id))
call neasyf_error(nf90_get_var(file_id, var_id, expected_data))
call neasyf_close(file_id)

@assertEqual(expected_data, read_data)
end subroutine test_write_scalar_unlimited

@test
subroutine test_write_1d
integer, parameter :: NX = 5
Expand Down

0 comments on commit 13f188d

Please sign in to comment.