submodule (penalisation_m) penalisation_netcdf_s !! Routines related with I/O of penalisation implicit none contains module subroutine write_netcdf_penalisation(self, fgid) class(penalisation_t), intent(in) :: self integer, intent(in) :: fgid integer :: n_points_grid integer :: nf90_stat, iddim_npoints_grid, id_charfun, id_dirindfun, & iddim_n_p_inds, iddim_n_pb_inds, id_p_inds, id_pb_inds ! define dimensions n_points_grid = size(self%charfun) nf90_stat = nf90_def_dim(fgid, 'n_points_grid', n_points_grid , iddim_npoints_grid) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_def_dim(fgid, 'n_p_inds', self%n_p_inds, iddim_n_p_inds) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_def_dim(fgid, 'n_pb_inds', self%n_pb_inds, iddim_n_pb_inds) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) ! define variables nf90_stat = nf90_def_var(fgid, 'charfun', NF90_DOUBLE, iddim_npoints_grid, id_charfun) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_def_var(fgid, 'dirindfun', NF90_DOUBLE, iddim_npoints_grid, id_dirindfun) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_def_var(fgid, 'p_inds', NF90_INT, iddim_n_p_inds, id_p_inds) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_def_var(fgid, 'pb_inds', NF90_INT, iddim_n_pb_inds, id_pb_inds) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) ! write global attributes (metadata) nf90_stat = nf90_put_att(fgid, NF90_GLOBAL, 'pen_method', self%pen_method) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_put_att(fgid, NF90_GLOBAL, 'penfuns_type', self%penfuns_type) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_put_att(fgid, NF90_GLOBAL, 'hermite_order', self%hermite_order) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_put_att(fgid, NF90_GLOBAL, 'charfun_parwidth', self%charfun_parwidth) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_put_att(fgid, NF90_GLOBAL, 'charfun_radlimwidth', self%charfun_radlimwidth) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_put_att(fgid, NF90_GLOBAL, 'dirindfun_parwidth', self%dirindfun_parwidth) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_put_att(fgid, NF90_GLOBAL, 'dphi', self%dphi) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_put_att(fgid, NF90_GLOBAL, 'dphi_max', self%dphi_max) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_put_att(fgid, NF90_GLOBAL, 'charfun_at_target', self%charfun_at_target) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) ! 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_charfun, self%charfun ) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_put_var(fgid, id_dirindfun, self%dirindfun ) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_put_var(fgid, id_p_inds, self%p_inds) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_put_var(fgid, id_pb_inds, self%pb_inds) 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_penalisation(self, fgid) class(penalisation_t), intent(inout) :: self integer, intent(in) :: fgid integer :: n_points_grid integer :: nf90_stat, iddim_npoints_grid, id_charfun, id_dirindfun, & iddim_n_p_inds, iddim_n_pb_inds, id_p_inds, id_pb_inds ! read global attributes nf90_stat = nf90_get_att(fgid, NF90_GLOBAL, 'pen_method', self%pen_method) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_get_att(fgid, NF90_GLOBAL, 'penfuns_type', self%penfuns_type) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_get_att(fgid, NF90_GLOBAL, 'hermite_order', self%hermite_order) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_get_att(fgid, NF90_GLOBAL, 'charfun_parwidth', self%charfun_parwidth) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_get_att(fgid, NF90_GLOBAL, 'charfun_radlimwidth', self%charfun_radlimwidth) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_get_att(fgid, NF90_GLOBAL, 'penfuns_type', self%penfuns_type) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_get_att(fgid, NF90_GLOBAL, 'dirindfun_parwidth', self%dirindfun_parwidth) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_get_att(fgid, NF90_GLOBAL, 'dphi', self%dphi) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_get_att(fgid, NF90_GLOBAL, 'dphi_max', self%dphi_max) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_get_att(fgid, NF90_GLOBAL, 'charfun_at_target', self%charfun_at_target) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) ! get dimensions nf90_stat = nf90_inq_dimid(fgid, 'n_points_grid', iddim_npoints_grid) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_inquire_dimension(fgid, iddim_npoints_grid, len = n_points_grid) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_inq_dimid(fgid, 'n_p_inds', iddim_n_p_inds) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_inquire_dimension(fgid, iddim_n_p_inds, len = self%n_p_inds) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_inq_dimid(fgid, 'n_pb_inds', iddim_n_pb_inds) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_inquire_dimension(fgid, iddim_n_pb_inds, len = self%n_pb_inds) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) ! allocate fields if (.not. allocated(self%charfun)) then allocate( self%charfun(n_points_grid) ) else call handle_error('Already allocated', & GRILLIX_ERR_ALLOC, __LINE__, __FILE__) endif if (.not. allocated(self%dirindfun)) then allocate( self%dirindfun(n_points_grid) ) else call handle_error('Already allocated', & GRILLIX_ERR_ALLOC, __LINE__, __FILE__) endif if ( .not. allocated(self%p_inds) ) then allocate(self%p_inds(self%n_p_inds)) else call handle_error('Already allocated', & GRILLIX_ERR_ALLOC, __LINE__, __FILE__) endif if ( .not. allocated(self%pb_inds) ) then allocate(self%pb_inds(self%n_pb_inds)) else call handle_error('Already allocated', & GRILLIX_ERR_ALLOC, __LINE__, __FILE__) endif ! read fields nf90_stat = nf90_inq_varid(fgid, 'charfun', id_charfun) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_get_var(fgid, id_charfun, self%charfun) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_inq_varid(fgid, 'dirindfun', id_dirindfun) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_get_var(fgid, id_dirindfun, self%dirindfun) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_inq_varid(fgid, 'p_inds', id_p_inds) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_get_var(fgid, id_p_inds, self%p_inds) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_inq_varid(fgid, 'pb_inds', id_pb_inds) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_get_var(fgid, id_pb_inds, self%pb_inds) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) end subroutine end submodule