submodule (parallel_target_flux_m) parallel_target_flux_netcdf_s !! Routines for parallel target marker file I/O implicit none contains module subroutine write_netcdf_target_markers(self, fgid) class(parallel_target_flux_t), intent(in) :: self integer, intent(in) :: fgid integer :: int_temp integer :: nf90_stat, iddim_nfwd, id_inds_fwd, id_weights_fwd, & iddim_nbwd, id_inds_bwd, id_weights_bwd ! Write metadata into global attributes int_temp = logical_to_integer(self%is_staggered) nf90_stat = nf90_put_att(fgid, NF90_GLOBAL, 'is_staggered', int_temp) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) ! Define dimensions nf90_stat = nf90_def_dim(fgid, 'n_inds_fwd', self%n_inds_fwd, & iddim_nfwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_def_dim(fgid, 'n_inds_bwd', self%n_inds_bwd, & iddim_nbwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) ! Define variables nf90_stat = nf90_def_var(fgid, 'inds_fwd', NF90_INT, iddim_nfwd, & id_inds_fwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_def_var(fgid, 'inds_bwd', NF90_INT, iddim_nbwd, & id_inds_bwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_def_var(fgid, 'weights_fwd', NF90_DOUBLE, iddim_nfwd, & id_weights_fwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_def_var(fgid, 'weights_bwd', NF90_DOUBLE, iddim_nbwd, & id_weights_bwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) ! End of definition phase nf90_stat = nf90_enddef(fgid) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) ! Put variables nf90_stat = nf90_put_var(fgid, id_inds_fwd, self%inds_fwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_put_var(fgid, id_inds_bwd, self%inds_bwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_put_var(fgid, id_weights_fwd, self%weights_fwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_put_var(fgid, id_weights_bwd, self%weights_bwd) 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_target_markers(self, fgid) class(parallel_target_flux_t), intent(inout) :: self integer, intent(in) :: fgid integer :: int_temp integer :: nf90_stat, iddim_nfwd, id_inds_fwd, id_weights_fwd, & iddim_nbwd, id_inds_bwd, id_weights_bwd ! Read metadata from global attributes nf90_stat = nf90_get_att(fgid, NF90_GLOBAL, 'is_staggered', int_temp) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) self%is_staggered = integer_to_logical(int_temp) ! Read dimensions nf90_stat = nf90_inq_dimid(fgid, 'n_inds_fwd', iddim_nfwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_inquire_dimension(fgid, iddim_nfwd, len = self%n_inds_fwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_inq_dimid(fgid, 'n_inds_bwd', iddim_nbwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_inquire_dimension(fgid, iddim_nbwd, len = self%n_inds_bwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) ! Allocate fields if ( .not. allocated(self%inds_fwd) ) then allocate(self%inds_fwd(self%n_inds_fwd)) else call handle_error('Already allocated', & GRILLIX_ERR_ALLOC, __LINE__, __FILE__) endif if ( .not. allocated(self%weights_fwd) ) then allocate(self%weights_fwd(self%n_inds_fwd)) else call handle_error('Already allocated', & GRILLIX_ERR_ALLOC, __LINE__, __FILE__) endif if ( .not. allocated(self%inds_bwd) ) then allocate(self%inds_bwd(self%n_inds_bwd)) else call handle_error('Already allocated', & GRILLIX_ERR_ALLOC, __LINE__, __FILE__) endif if ( .not. allocated(self%weights_bwd) ) then allocate(self%weights_bwd(self%n_inds_bwd)) else call handle_error('Already allocated', & GRILLIX_ERR_ALLOC, __LINE__, __FILE__) endif ! Read fields nf90_stat = nf90_inq_varid(fgid, 'inds_fwd', id_inds_fwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_get_var(fgid, id_inds_fwd, self%inds_fwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_inq_varid(fgid, 'weights_fwd', id_weights_fwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_get_var(fgid, id_weights_fwd, self%weights_fwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_inq_varid(fgid, 'inds_bwd', id_inds_bwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_get_var(fgid, id_inds_bwd, self%inds_bwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_inq_varid(fgid, 'weights_bwd', id_weights_bwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) nf90_stat = nf90_get_var(fgid, id_weights_bwd, self%weights_bwd) call handle_error_netcdf(nf90_stat, __LINE__, __FILE__) end subroutine end submodule