parbnd_taylor_netcdf_s.f90 Source File


This file depends on

sourcefile~~parbnd_taylor_netcdf_s.f90~~EfferentGraph sourcefile~parbnd_taylor_netcdf_s.f90 parbnd_taylor_netcdf_s.f90 sourcefile~error_handling_m.f90 error_handling_m.f90 sourcefile~parbnd_taylor_netcdf_s.f90->sourcefile~error_handling_m.f90 sourcefile~parbnd_taylor_m.f90 parbnd_taylor_m.f90 sourcefile~parbnd_taylor_netcdf_s.f90->sourcefile~parbnd_taylor_m.f90 sourcefile~comm_handling_m.f90 comm_handling_m.f90 sourcefile~error_handling_m.f90->sourcefile~comm_handling_m.f90 sourcefile~precision_m.f90 precision_m.f90 sourcefile~error_handling_m.f90->sourcefile~precision_m.f90 sourcefile~screen_io_m.f90 screen_io_m.f90 sourcefile~error_handling_m.f90->sourcefile~screen_io_m.f90 sourcefile~status_codes_m.f90 status_codes_m.f90 sourcefile~error_handling_m.f90->sourcefile~status_codes_m.f90 sourcefile~parbnd_taylor_m.f90->sourcefile~error_handling_m.f90 sourcefile~connection_length_m.f90 connection_length_m.f90 sourcefile~parbnd_taylor_m.f90->sourcefile~connection_length_m.f90 sourcefile~constants_m.f90 constants_m.f90 sourcefile~parbnd_taylor_m.f90->sourcefile~constants_m.f90 sourcefile~elementary_functions_m.f90 elementary_functions_m.f90 sourcefile~parbnd_taylor_m.f90->sourcefile~elementary_functions_m.f90 sourcefile~equilibrium_m.f90 equilibrium_m.f90 sourcefile~parbnd_taylor_m.f90->sourcefile~equilibrium_m.f90 sourcefile~fieldline_tracer_m.f90 fieldline_tracer_m.f90 sourcefile~parbnd_taylor_m.f90->sourcefile~fieldline_tracer_m.f90 sourcefile~mesh_cart_m.f90 mesh_cart_m.f90 sourcefile~parbnd_taylor_m.f90->sourcefile~mesh_cart_m.f90 sourcefile~parbnd_taylor_m.f90->sourcefile~precision_m.f90 sourcefile~parbnd_taylor_m.f90->sourcefile~screen_io_m.f90 sourcefile~parbnd_taylor_m.f90->sourcefile~status_codes_m.f90 sourcefile~connection_length_m.f90->sourcefile~error_handling_m.f90 sourcefile~connection_length_m.f90->sourcefile~equilibrium_m.f90 sourcefile~connection_length_m.f90->sourcefile~fieldline_tracer_m.f90 sourcefile~connection_length_m.f90->sourcefile~precision_m.f90 sourcefile~connection_length_m.f90->sourcefile~screen_io_m.f90 sourcefile~connection_length_m.f90->sourcefile~status_codes_m.f90 sourcefile~descriptors_m.f90 descriptors_m.f90 sourcefile~connection_length_m.f90->sourcefile~descriptors_m.f90 sourcefile~constants_m.f90->sourcefile~precision_m.f90 sourcefile~elementary_functions_m.f90->sourcefile~precision_m.f90 sourcefile~equilibrium_m.f90->sourcefile~precision_m.f90 sourcefile~fieldline_tracer_m.f90->sourcefile~error_handling_m.f90 sourcefile~fieldline_tracer_m.f90->sourcefile~comm_handling_m.f90 sourcefile~fieldline_tracer_m.f90->sourcefile~equilibrium_m.f90 sourcefile~fieldline_tracer_m.f90->sourcefile~precision_m.f90 sourcefile~fieldline_tracer_m.f90->sourcefile~screen_io_m.f90 sourcefile~fieldline_tracer_m.f90->sourcefile~status_codes_m.f90 sourcefile~mesh_cart_m.f90->sourcefile~error_handling_m.f90 sourcefile~mesh_cart_m.f90->sourcefile~comm_handling_m.f90 sourcefile~mesh_cart_m.f90->sourcefile~equilibrium_m.f90 sourcefile~mesh_cart_m.f90->sourcefile~precision_m.f90 sourcefile~mesh_cart_m.f90->sourcefile~status_codes_m.f90 sourcefile~mesh_cart_m.f90->sourcefile~descriptors_m.f90 sourcefile~slab_equilibrium.f90 slab_equilibrium.f90 sourcefile~mesh_cart_m.f90->sourcefile~slab_equilibrium.f90 sourcefile~screen_io_m.f90->sourcefile~precision_m.f90 sourcefile~descriptors_m.f90->sourcefile~screen_io_m.f90 sourcefile~slab_equilibrium.f90->sourcefile~equilibrium_m.f90 sourcefile~slab_equilibrium.f90->sourcefile~precision_m.f90 sourcefile~slab_equilibrium.f90->sourcefile~screen_io_m.f90 sourcefile~slab_equilibrium.f90->sourcefile~descriptors_m.f90 sourcefile~params_equi_slab_m.f90 params_equi_slab_m.f90 sourcefile~slab_equilibrium.f90->sourcefile~params_equi_slab_m.f90 sourcefile~params_equi_slab_m.f90->sourcefile~error_handling_m.f90 sourcefile~params_equi_slab_m.f90->sourcefile~precision_m.f90 sourcefile~params_equi_slab_m.f90->sourcefile~screen_io_m.f90 sourcefile~params_equi_slab_m.f90->sourcefile~status_codes_m.f90

Source Code

submodule(parbnd_taylor_m) parbnd_taylor_netcdf_s
    !! Routines related with NetCDF I/O of parbnd_taylor
    use netcdf
    use error_handling_m, only : handle_error_netcdf
    implicit none

contains

    module subroutine write_netcdf_parbnd_taylor(self, fgid)
        class(parbnd_taylor_t), intent(in) :: self
        integer, intent(in) :: fgid

        integer :: np_mesh, nf90_stat
        integer :: iddim_np_mesh, iddim_max_taylor_order, &
            id_b_bnd_direction, id_dists_to_bnd, &
            id_depth_in_bnd

        ! Define dimensions 
        np_mesh = size(self%b_bnd_direction)
        nf90_stat = nf90_def_dim(fgid, 'np_mesh', np_mesh , iddim_np_mesh)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_def_dim(fgid, 'max_taylor_order', &
                max_taylor_order, iddim_max_taylor_order)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        ! Define variables
        nf90_stat = nf90_def_var(fgid, 'b_bnd_direction', NF90_INT, &
            iddim_np_mesh, id_b_bnd_direction)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_def_var(fgid, 'dists_to_bnd', NF90_DOUBLE, &
            (/ iddim_np_mesh, iddim_max_taylor_order/), id_dists_to_bnd)

        nf90_stat = nf90_def_var(fgid, 'depth_in_bnd', NF90_DOUBLE, &
            (/ iddim_np_mesh /), id_depth_in_bnd)

        ! End of definition
        nf90_stat = nf90_enddef(fgid)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        ! Put variables
        nf90_stat = nf90_put_var(fgid, id_b_bnd_direction, &
            self%b_bnd_direction)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_put_var(fgid, id_dists_to_bnd, &
            self%dists_to_bnd)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_put_var(fgid, id_depth_in_bnd, &
            self%depth_in_bnd)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        ! put into redefinition state
        nf90_stat = nf90_redef(fgid)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

    end subroutine

    module subroutine read_netcdf_parbnd_taylor(self, fgid)
        class(parbnd_taylor_t), intent(inout) :: self
        integer, intent(in) :: fgid

        integer :: np_mesh, nf90_stat
        integer :: iddim_np_mesh, iddim_max_taylor_order, &
            id_b_bnd_direction, id_dists_to_bnd, &
            id_depth_in_bnd
        integer :: max_taylor_order_in

        ! Get dimensions
        nf90_stat = nf90_inq_dimid(fgid, 'np_mesh', iddim_np_mesh)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        nf90_stat =  nf90_inquire_dimension(fgid, iddim_np_mesh, len=np_mesh)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_inq_dimid(fgid, 'max_taylor_order', &
            iddim_max_taylor_order)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        nf90_stat =  nf90_inquire_dimension(fgid, iddim_max_taylor_order, &
            len=max_taylor_order_in)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        ! Check if order max_taylor_order from file agrees
        if (max_taylor_order_in /= max_taylor_order) then
            call handle_error('Mismatch of max_taylor_order between netcdf &
                               file and code', PARALLAX_ERR_PARBND, &
                              __LINE__, __FILE__)
        endif

        ! Allocate fields
        allocate(self%b_bnd_direction(np_mesh))
        allocate(self%dists_to_bnd(np_mesh, max_taylor_order))
        allocate(self%depth_in_bnd(np_mesh))

        ! Read fields
        nf90_stat = nf90_inq_varid(fgid, 'b_bnd_direction', &
            id_b_bnd_direction)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        nf90_stat = nf90_get_var(fgid, id_b_bnd_direction, &
            self%b_bnd_direction)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_inq_varid(fgid, 'dists_to_bnd', &
            id_dists_to_bnd)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        nf90_stat = nf90_get_var(fgid, id_dists_to_bnd, &
            self%dists_to_bnd)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

        nf90_stat = nf90_inq_varid(fgid, 'depth_in_bnd', &
            id_depth_in_bnd)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)
        nf90_stat = nf90_get_var(fgid, id_depth_in_bnd, &
            self%depth_in_bnd)
        call handle_error_netcdf(nf90_stat, __LINE__, __FILE__)

    end subroutine

end submodule