module netcdf_io_module
    use netcdf
    implicit none

    type :: io_struct
        integer :: ncid
        integer :: time_dimid, height_dimid, lat_dimid, lon_dimid
        integer :: varid_time, varid_series, varid_timescan, varid_profile
        integer :: varid_2Dtime, varid_3Dtime, varid_2D, varid_3D
        logical :: is_open = .false.
        character(len=128) :: series_name = 'timeseries of variable'
        character(len=128) :: series_long_name = 'timeseries of some variable'
        character(len=128) :: series_units = 'unit'
        character(len=128) :: timescan_name = 'vertical timescan of variable'
        character(len=128) :: timescan_long_name = 'vertical timescan of some variable'
        character(len=128) :: timescan_units = 'unit'
        character(len=128) :: profile_name = 'vertical profile of variable'
        character(len=128) :: profile_long_name = 'vertical profile of some variable'
        character(len=128) :: profile_units = 'unit'
        character(len=128) :: d2time_name = '2D variable with time'
        character(len=128) :: d2time_long_name = '2D variable with time'
        character(len=128) :: d2time_units = 'unit'
        character(len=128) :: d3time_name = '3D variable with time'
        character(len=128) :: d3time_long_name = '3D variable with time'
        character(len=128) :: d3time_units = 'unit'
        character(len=128) :: d2_name = '2D variable'
        character(len=128) :: d2_long_name = '2D variable'
        character(len=128) :: d2_units = 'unit'
        character(len=128) :: d3_name = '3D variable'
        character(len=128) :: d3_long_name = '3D variable'
        character(len=128) :: d3_units = 'unit'
    end type io_struct

contains
    subroutine open_netcdf(filename, ios, time_len, height_len, lat_len, lon_len)
        character(len=*), intent(in) :: filename
        type(io_struct), intent(out) :: ios
        integer, intent(in) :: time_len, height_len, lat_len, lon_len
        integer :: ierr

        ierr = nf90_create(filename, nf90_clobber, ios%ncid)
        if (ierr == nf90_noerr) then
            ios%is_open = .true.
            ierr = nf90_def_dim(ios%ncid, 'time', time_len, ios%time_dimid)
            ierr = nf90_def_dim(ios%ncid, 'height', height_len, ios%height_dimid)
            ierr = nf90_def_dim(ios%ncid, 'lat', lat_len, ios%lat_dimid)
            ierr = nf90_def_dim(ios%ncid, 'lon', lon_len, ios%lon_dimid)
            ierr = nf90_def_var(ios%ncid, 'time', nf90_float, ios%time_dimid, ios%varid_time)
            ierr = nf90_def_var(ios%ncid, 'series', nf90_float, ios%time_dimid, ios%varid_series)
            ierr = nf90_def_var(ios%ncid, 'timescan', nf90_float, &
                    [ios%time_dimid, ios%height_dimid], ios%varid_timescan)
            ierr = nf90_def_var(ios%ncid, 'profile', nf90_float, ios%height_dimid, ios%varid_profile)
            ierr = nf90_def_var(ios%ncid, '2Dtime', nf90_float, &
                    [ios%time_dimid, ios%lat_dimid, ios%lon_dimid], ios%varid_2Dtime)
            ierr = nf90_def_var(ios%ncid, '3Dtime', nf90_float, &
                    [ios%time_dimid, ios%lat_dimid, ios%lon_dimid, ios%height_dimid], ios%varid_3Dtime)
            ierr = nf90_def_var(ios%ncid, '2D', nf90_float, &
                    [ios%lat_dimid, ios%lon_dimid], ios%varid_2D)
            ierr = nf90_def_var(ios%ncid, '3D', nf90_float, &
                    [ios%lat_dimid, ios%lon_dimid, ios%height_dimid], ios%varid_3D)

            ierr = nf90_put_att(ios%ncid, ios%varid_series, 'standard_name', ios%series_name)
            ierr = nf90_put_att(ios%ncid, ios%varid_series, 'long_name', ios%series_long_name)
            ierr = nf90_put_att(ios%ncid, ios%varid_series, 'units', ios%series_units)
            ierr = nf90_put_att(ios%ncid, ios%varid_timescan, 'standard_name', ios%timescan_name)
            ierr = nf90_put_att(ios%ncid, ios%varid_timescan, 'long_name', ios%timescan_long_name)
            ierr = nf90_put_att(ios%ncid, ios%varid_timescan, 'units', ios%timescan_units)
            ierr = nf90_put_att(ios%ncid, ios%varid_profile, 'standard_name', ios%profile_name)
            ierr = nf90_put_att(ios%ncid, ios%varid_profile, 'long_name', ios%profile_long_name)
            ierr = nf90_put_att(ios%ncid, ios%varid_profile, 'units', ios%profile_units)
            ierr = nf90_put_att(ios%ncid, ios%varid_2Dtime, 'standard_name', ios%d2time_name)
            ierr = nf90_put_att(ios%ncid, ios%varid_2Dtime, 'long_name', ios%d2time_long_name)
            ierr = nf90_put_att(ios%ncid, ios%varid_2Dtime, 'units', ios%d2time_units)
            ierr = nf90_put_att(ios%ncid, ios%varid_3Dtime, 'standard_name', ios%d3time_name)
            ierr = nf90_put_att(ios%ncid, ios%varid_3Dtime, 'long_name', ios%d3time_long_name)
            ierr = nf90_put_att(ios%ncid, ios%varid_3Dtime, 'units', ios%d3time_units)
            ierr = nf90_put_att(ios%ncid, ios%varid_2D, 'standard_name', ios%d2_name)
            ierr = nf90_put_att(ios%ncid, ios%varid_2D, 'long_name', ios%d2_long_name)
            ierr = nf90_put_att(ios%ncid, ios%varid_2D, 'units', ios%d2_units)
            ierr = nf90_put_att(ios%ncid, ios%varid_3D, 'standard_name', ios%d3_name)
            ierr = nf90_put_att(ios%ncid, ios%varid_3D, 'long_name', ios%d3_long_name)
            ierr = nf90_put_att(ios%ncid, ios%varid_3D, 'units', ios%d3_units)

            ierr = nf90_enddef(ios%ncid)
        else
            print *, 'Error opening file: ', ierr
            ios%is_open = .false.
        endif
    end subroutine open_netcdf

    subroutine write_series(ios, time_data, series_data)
        type(io_struct), intent(in) :: ios
        real, dimension(:), intent(in) :: time_data, series_data
        integer :: ierr

        if (.not. ios%is_open) then
            print *, "File is not open."
            return
        endif

        ierr = nf90_put_var(ios%ncid, ios%varid_time, time_data)
        ierr = nf90_put_var(ios%ncid, ios%varid_series, series_data)
        if (ierr /= nf90_noerr) then
            print *, 'Error writing series data:', ierr
        endif
    end subroutine write_series

    subroutine write_timescan(ios, timescan_data)
        type(io_struct), intent(in) :: ios
        real, dimension(:,:), intent(in) :: timescan_data
        integer :: ierr

        if (.not. ios%is_open) then
            print *, "File is not open."
            return
        endif

        ierr = nf90_put_var(ios%ncid, ios%varid_timescan, timescan_data)
        if (ierr /= nf90_noerr) then
            print *, 'Error writing timescan data:', ierr
        endif
    end subroutine write_timescan

    subroutine write_profile(ios, profile_data)
        type(io_struct), intent(in) :: ios
        real, dimension(:), intent(in) :: profile_data
        integer :: ierr

        if (.not. ios%is_open) then
            print *, "File is not open."
            return
        endif

        ierr = nf90_put_var(ios%ncid, ios%varid_profile, profile_data)
        if (ierr /= nf90_noerr) then
            print *, 'Error writing profile data:', ierr
        endif
    end subroutine write_profile

    subroutine write_2Dtime(ios, d2time_data)
        type(io_struct), intent(in) :: ios
        real, dimension(:,:,:) :: d2time_data
        integer :: ierr

        if (.not. ios%is_open) then
            print *, "File is not open."
            return
        endif

        ierr = nf90_put_var(ios%ncid, ios%varid_2Dtime, d2time_data)
        if (ierr /= nf90_noerr) then
            print *, 'Error writing 2Dtime data:', ierr
        endif
    end subroutine write_2Dtime

    subroutine write_3Dtime(ios, d3time_data)
        type(io_struct), intent(in) :: ios
        real, dimension(:,:,:,:) :: d3time_data
        integer :: ierr

        if (.not. ios%is_open) then
            print *, "File is not open."
            return
        endif

        ierr = nf90_put_var(ios%ncid, ios%varid_3Dtime, d3time_data)
        if (ierr /= nf90_noerr) then
            print *, "Error writing 3Dtime data:", ierr
        endif
    end subroutine write_3Dtime

    subroutine write_2D(ios, d2_data)
        type(io_struct), intent(in) :: ios
        real, dimension(:,:) :: d2_data
        integer :: ierr

        if (.not. ios%is_open) then
            print *, "File is not open."
            return
        endif

        ierr = nf90_put_var(ios%ncid, ios%varid_2D, d2_data)
        if (ierr /= nf90_noerr) then
            print *, 'Error writing 2D data:', ierr
        endif
    end subroutine write_2D

    subroutine write_3D(ios, d3_data)
        type(io_struct), intent(in) :: ios
        real, dimension(:,:,:) :: d3_data
        integer :: ierr

        if (.not. ios%is_open) then
            print *, "File is not open."
            return
        endif

        ierr = nf90_put_var(ios%ncid, ios%varid_3D, d3_data)
        if (ierr /= nf90_noerr) then
            print *, 'Error writing 3D data:', ierr
        endif
    end subroutine write_3D

    subroutine close_netcdf(ios)
        type(io_struct), intent(inout) :: ios
        integer :: ierr

        ierr = nf90_close(ios%ncid)
        ios%is_open = .false.
    end subroutine close_netcdf
end module netcdf_io_module