Skip to content

Commit

Permalink
add and check is_open property
Browse files Browse the repository at this point in the history
  • Loading branch information
scivision committed Jul 16, 2020
1 parent 0ab8d66 commit 3bb648d
Show file tree
Hide file tree
Showing 5 changed files with 48 additions and 2 deletions.
7 changes: 6 additions & 1 deletion src/interface.f90
Original file line number Diff line number Diff line change
Expand Up @@ -224,7 +224,7 @@ subroutine nc_initialize(self,filename,ierr, status,action,comp_lvl,verbose,debu
integer :: ier

if (self%is_open) then
write(stderr,*) 'WARNING:nc4fortran: file handle already open to: '// filename
write(stderr,*) 'WARNING:nc4fortran:initialize file handle already open to: '// filename
return
endif

Expand Down Expand Up @@ -287,6 +287,11 @@ subroutine nc_finalize(self, ierr)

integer :: ier

if(.not. self%is_open) then
write(stderr,*) 'WARNING:nc4fortran:finalize file handle is not open'
return
endif

ier = nf90_close(self%ncid)
if (present(ierr)) ierr = ier
if (ier /= NF90_NOERR) then
Expand Down
13 changes: 12 additions & 1 deletion src/read.f90
Original file line number Diff line number Diff line change
Expand Up @@ -11,6 +11,9 @@
integer :: i, tempdims(NC_MAXDIM), N, ier
character(NF90_MAX_NAME) :: tempnames(NC_MAXDIM)


if(.not.self%is_open) error stop 'ERROR:nc4fortran:read: file handle not open'

N = 0
do i = 1,NC_MAXDIM
ier = nf90_inquire_dimension(self%ncid, dimid=i, name=tempnames(i), len=tempdims(i))
Expand All @@ -28,9 +31,17 @@

module procedure nc_check_exist
integer :: varid, ierr
ierr = nf90_inq_varid(self%ncid, dname, varid)

exists = .false.

if(.not.self%is_open) then
write(stderr,*) 'WARNING:nc4fortran:read: file handle not open ' // self%filename
return
endif

ierr = nf90_inq_varid(self%ncid, dname, varid)


select case (ierr)
case (NF90_NOERR)
exists = .true.
Expand Down
24 changes: 24 additions & 0 deletions src/reader.f90
Original file line number Diff line number Diff line change
Expand Up @@ -7,6 +7,9 @@

module procedure nc_read_scalar
integer :: varid, ier

if(.not.self%is_open) error stop 'ERROR:nc4fortran:reader file handle not open'

ier = nf90_inq_varid(self%ncid, dname, varid)

if(ier == NF90_NOERR) then
Expand Down Expand Up @@ -36,6 +39,9 @@

module procedure nc_read_1d
integer :: varid, ier

if(.not.self%is_open) error stop 'ERROR:nc4fortran:reader file handle not open'

ier = nf90_inq_varid(self%ncid, dname, varid)

if(ier == NF90_NOERR) then
Expand Down Expand Up @@ -65,6 +71,9 @@

module procedure nc_read_2d
integer :: varid, ier

if(.not.self%is_open) error stop 'ERROR:nc4fortran:reader file handle not open'

ier = nf90_inq_varid(self%ncid, dname, varid)

if(ier == NF90_NOERR) then
Expand Down Expand Up @@ -93,6 +102,9 @@

module procedure nc_read_3d
integer :: varid, ier

if(.not.self%is_open) error stop 'ERROR:nc4fortran:reader file handle not open'

ier = nf90_inq_varid(self%ncid, dname, varid)

if(ier == NF90_NOERR) then
Expand Down Expand Up @@ -121,6 +133,9 @@

module procedure nc_read_4d
integer :: varid, ier

if(.not.self%is_open) error stop 'ERROR:nc4fortran:reader file handle not open'

ier = nf90_inq_varid(self%ncid, dname, varid)

if(ier == NF90_NOERR) then
Expand Down Expand Up @@ -149,6 +164,9 @@

module procedure nc_read_5d
integer :: varid, ier

if(.not.self%is_open) error stop 'ERROR:nc4fortran:reader file handle not open'

ier = nf90_inq_varid(self%ncid, dname, varid)

if(ier == NF90_NOERR) then
Expand Down Expand Up @@ -177,6 +195,9 @@

module procedure nc_read_6d
integer :: varid, ier

if(.not.self%is_open) error stop 'ERROR:nc4fortran:reader file handle not open'

ier = nf90_inq_varid(self%ncid, dname, varid)

if(ier == NF90_NOERR) then
Expand Down Expand Up @@ -205,6 +226,9 @@

module procedure nc_read_7d
integer :: varid, ier

if(.not.self%is_open) error stop 'ERROR:nc4fortran:reader file handle not open'

ier = nf90_inq_varid(self%ncid, dname, varid)

if(ier == NF90_NOERR) then
Expand Down
4 changes: 4 additions & 0 deletions src/write.f90
Original file line number Diff line number Diff line change
Expand Up @@ -10,6 +10,8 @@
integer :: i, L
character(NF90_MAX_NAME) :: name

if(.not.self%is_open) error stop 'ERROR:nc4fortran:write: file handle not open'

do i=1,size(dims)
if (present(dimnames)) then
ierr = nf90_inq_dimid(self%ncid, dimnames(i), dimids(i))
Expand All @@ -34,6 +36,8 @@
module procedure write_attribute
integer :: varid, ier

if(.not.self%is_open) error stop 'ERROR:nc4fortran:write: file handle not open'

ier = nf90_inq_varid(self%ncid, dname, varid)

if(ierr == nf90_noerr) ier = nf90_put_att(self%ncid, varid, attrname, value)
Expand Down
2 changes: 2 additions & 0 deletions src/writer.f90
Original file line number Diff line number Diff line change
Expand Up @@ -8,6 +8,8 @@
module procedure nc_write_scalar
integer :: varid, ier

if(.not.self%is_open) error stop 'ERROR:nc4fortran:writer: file handle not open'

select type (value)
type is (real(real64))
ier = nf90_def_var(self%ncid, dname, NF90_DOUBLE, varid=varid)
Expand Down

0 comments on commit 3bb648d

Please sign in to comment.