Skip to content
Snippets Groups Projects
Commit e14c1952 authored by GUEZ Lionel's avatar GUEZ Lionel
Browse files

Move procedure `get_var` to its own module

parent e8d1d6b9
No related branches found
No related tags found
No related merge requests found
# inst_eddies:
add_executable(inst_eddies inst_eddies.f90 local_extrema.f90 set_max_speed.f90
get_1_outerm.f90 good_contour.f90 mean_speed.f90 set_contours.f90
nearby_extr.f90 config.f90 input_ssh.f90 read_grid_inst_eddies.f90
nearby_extr.f90 config.f90 input_ssh.f90 get_var.f90 read_grid_inst_eddies.f90
write_snapshot.f90 ccw_orient.f90 write_eddy.f90 set_all_extr.f90
complete_ssh.f90 invert_proj.f90)
target_link_libraries(inst_eddies PRIVATE Contour_531::contour_531
......@@ -26,18 +26,18 @@ configure_file(find_coord.py . COPYONLY)
# Tests:
add_subdirectory(Tests)
target_sources(test_get_1_outerm PRIVATE get_1_outerm.f90 good_contour.f90
config.f90 input_ssh.f90 read_grid_inst_eddies.f90 write_eddy.f90
config.f90 input_ssh.f90 get_var.f90 read_grid_inst_eddies.f90 write_eddy.f90
ccw_orient.f90 invert_proj.f90)
target_sources(test_good_contour PRIVATE good_contour.f90)
target_sources(test_mean_speed PRIVATE mean_speed.f90 input_ssh.f90
target_sources(test_mean_speed PRIVATE mean_speed.f90 input_ssh.f90 get_var.f90
read_grid_inst_eddies.f90 config.f90)
target_sources(test_nearby_extr PRIVATE nearby_extr.f90)
target_sources(test_local_extrema PRIVATE local_extrema.f90)
target_sources(test_set_max_speed PRIVATE set_max_speed.f90 good_contour.f90
mean_speed.f90 config.f90 input_ssh.f90 read_grid_inst_eddies.f90
mean_speed.f90 config.f90 input_ssh.f90 get_var.f90 read_grid_inst_eddies.f90
ccw_orient.f90 write_eddy.f90 complete_ssh.f90 invert_proj.f90)
target_sources(test_complete_ssh PRIVATE complete_ssh.f90)
target_sources(examine_eddy PRIVATE config.f90 input_ssh.f90
target_sources(examine_eddy PRIVATE config.f90 input_ssh.f90 get_var.f90
read_grid_inst_eddies.f90 set_all_extr.f90 local_extrema.f90 write_eddy.f90
set_contours.f90 get_1_outerm.f90 ccw_orient.f90 good_contour.f90
set_max_speed.f90 complete_ssh.f90 mean_speed.f90 invert_proj.f90)
module get_var_m
implicit none
contains
subroutine get_var(periodic, values, ncid, nc_time, name, new_fill_value)
! Read a NetCDF variable, change the missing value and extend it
! in longitude if periodic.
! This subroutine is in module input_ssh_m because it needs to
! access nlon, and it is used by input_ssh, so it cannot use the
! module input_ssh_m, which would be a circular dependency.
! Libraries:
use netcdf95, only: nf95_inq_varid, nf95_gunp_var
use config_m, only: max_radius
use read_grid_inst_eddies_m, only: nlon
logical, intent(in):: periodic ! grid is periodic in longitude
real, intent(out):: values(1 - merge(max_radius(1), 0, periodic):, :)
! (1 - merge(max_radius(1), 0, periodic):nlon + merge(max_radius(1),
! 0, periodic), nlat)
! ssh, u or v. We cannot place this argument first because the
! declaration references periodic.
integer, intent(in):: ncid
integer, intent(in):: nc_time
! Index, 1-based, in the NetCDF time coordinate, if any, of the
! time slice to read. If there is no time coordinate, then nc_time
! should be equal to 1. If there is a time coordinate then we
! assume that the dimensions of the NetCDF variable to read are
! longitude, latitude and time, in Fortran order.
character(len = *), intent(in):: name ! of NetCDF variable
real, intent(in):: new_fill_value
! Local:
integer varid
!-------------------------------------------------------------------------
call nf95_inq_varid(ncid, name, varid)
if (nc_time == 1) then
call nf95_gunp_var(ncid, varid, values(1:nlon, :), &
new_missing = new_fill_value)
else
call nf95_gunp_var(ncid, varid, values(1:nlon, :), &
start = [1, 1, nc_time], new_missing = new_fill_value)
end if
if (periodic) then
! Extend in longitude:
values(1 - max_radius(1):0, :) &
= values(nlon + 1 - max_radius(1):nlon, :)
values(nlon + 1:nlon + max_radius(1), :) = values(1:max_radius(1), :)
end if
end subroutine get_var
end module get_var_m
......@@ -12,7 +12,6 @@ module input_ssh_m
! (1 - max_radius(1):nlon + max_radius(1), nlat) if periodic, else
! (nlon, nlat). Wind, in m s-1.
private get_var
contains
......@@ -25,6 +24,7 @@ contains
use netcdf95, only: nf95_open, nf95_close, nf95_nowrite
use config_m, only: max_radius
use get_var_m, only: get_var
use read_grid_inst_eddies_m, only: read_grid_inst_eddies, nlon, nlat
logical, intent(out):: periodic ! grid is periodic in longitude
......@@ -102,65 +102,4 @@ contains
end subroutine input_ssh
!*********************************************************************
subroutine get_var(periodic, values, ncid, nc_time, name, new_fill_value)
! Read a NetCDF variable, change the missing value and extend it
! in longitude if periodic.
! This subroutine is in module input_ssh_m because it needs to
! access nlon, and it is used by input_ssh, so it cannot use the
! module input_ssh_m, which would be a circular dependency.
! Libraries:
use netcdf95, only: nf95_inq_varid, nf95_gunp_var
use config_m, only: max_radius
use read_grid_inst_eddies_m, only: nlon
logical, intent(in):: periodic ! grid is periodic in longitude
real, intent(out):: values(1 - merge(max_radius(1), 0, periodic):, :)
! (1 - merge(max_radius(1), 0, periodic):nlon + merge(max_radius(1),
! 0, periodic), nlat)
! ssh, u or v. We cannot place this argument first because the
! declaration references periodic.
integer, intent(in):: ncid
integer, intent(in):: nc_time
! Index, 1-based, in the NetCDF time coordinate, if any, of the
! time slice to read. If there is no time coordinate, then nc_time
! should be equal to 1. If there is a time coordinate then we
! assume that the dimensions of the NetCDF variable to read are
! longitude, latitude and time, in Fortran order.
character(len = *), intent(in):: name ! of NetCDF variable
real, intent(in):: new_fill_value
! Local:
integer varid
!-------------------------------------------------------------------------
call nf95_inq_varid(ncid, name, varid)
if (nc_time == 1) then
call nf95_gunp_var(ncid, varid, values(1:nlon, :), &
new_missing = new_fill_value)
else
call nf95_gunp_var(ncid, varid, values(1:nlon, :), &
start = [1, 1, nc_time], new_missing = new_fill_value)
end if
if (periodic) then
! Extend in longitude:
values(1 - max_radius(1):0, :) &
= values(nlon + 1 - max_radius(1):nlon, :)
values(nlon + 1:nlon + max_radius(1), :) = values(1:max_radius(1), :)
end if
end subroutine get_var
end module input_ssh_m
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment