-
Lionel GUEZ authored
Pass argument copy to `read_snapshot` rather than arguments periodic and `dist_lim`. Motivation: do not recompute copy at each call of `read_snapshot`. `read_snapshot` is going to be called many times in overlap program, with always the same value of copy. And we do not need the values of periodic and `dist_lim` separately in `read_snapshot`. Complete `send_snapshot` and `recv_snapshot` by sending and receiving an entire snapshot. So we need nlon, nlat and copy in `recv_snapshot` to allocate `extr_map`. It is better to pass arguments nlon, nlat and copy than sending and receiving the shape of `extr_map` with MPI, of course, since the shape is the same for all snapshots. In program `test_send_recv`, write the entire received snaphot. Bug fix in `Tests/CMakeLists.txt`: several targets missed `${netcdff_INCLUDE_DIR}`. Create procedure `write_snapshot` from code in main program `test_read_snapshot`, because we need to do this also in program `test_send_recv`. Note the dissymmetry with `read_snapshot`, which does not include the opening and closing of the shapefiles. That is because we anticipate calling `read_snapshot` several times with the same shapefiles in program overlap. Change order of dummy arguments of procedure `write_extr_map`. More convenient for calling it with keywords from `write_snapshot`.
Lionel GUEZ authoredPass argument copy to `read_snapshot` rather than arguments periodic and `dist_lim`. Motivation: do not recompute copy at each call of `read_snapshot`. `read_snapshot` is going to be called many times in overlap program, with always the same value of copy. And we do not need the values of periodic and `dist_lim` separately in `read_snapshot`. Complete `send_snapshot` and `recv_snapshot` by sending and receiving an entire snapshot. So we need nlon, nlat and copy in `recv_snapshot` to allocate `extr_map`. It is better to pass arguments nlon, nlat and copy than sending and receiving the shape of `extr_map` with MPI, of course, since the shape is the same for all snapshots. In program `test_send_recv`, write the entire received snaphot. Bug fix in `Tests/CMakeLists.txt`: several targets missed `${netcdff_INCLUDE_DIR}`. Create procedure `write_snapshot` from code in main program `test_read_snapshot`, because we need to do this also in program `test_send_recv`. Note the dissymmetry with `read_snapshot`, which does not include the opening and closing of the shapefiles. That is because we anticipate calling `read_snapshot` several times with the same shapefiles in program overlap. Change order of dummy arguments of procedure `write_extr_map`. More convenient for calling it with keywords from `write_snapshot`.
write_extr_map.f90 1.60 KiB
module write_extr_map_m
implicit none
contains
subroutine write_extr_map(extr_map, longitude, latitude)
use netcdf, only: nf90_clobber, NF90_FLOAT, NF90_INT
use netcdf95, only: nf95_close, nf95_put_var, nf95_create, nf95_def_dim, &
nf95_def_var, nf95_put_att, nf95_enddef
integer, intent(in):: extr_map(:, :) ! map of extrema
real, intent(in):: longitude(:) ! in degrees
real, intent(in):: latitude(:) ! in degrees
! Local:
integer ncid, varid_lat, varid_lon, varid_extr_map, dimid_lat, dimid_lon
!--------------------------------------------------------------------------
call nf95_create("extr_map.nc", NF90_CLOBBER, ncid)
call nf95_def_dim(ncid, "lat", size(latitude), dimid_lat)
call nf95_def_dim(ncid, "lon", size(longitude), dimid_lon)
call nf95_def_var(ncid, "lat", NF90_FLOAT, dimid_lat, varid_lat)
call nf95_put_att(ncid, varid_lat, "standard_name", "latitude")
call nf95_put_att(ncid, varid_lat, "units", "degrees_north")
call nf95_def_var(ncid, "lon", NF90_FLOAT, dimid_lon, varid_lon)
call nf95_put_att(ncid, varid_lon, "standard_name", "longitude")
call nf95_put_att(ncid, varid_lon, "units", "degrees_east")
call nf95_def_var(ncid, "extr_map", NF90_INT, [dimid_lon, dimid_lat], &
varid_extr_map)
call nf95_enddef(ncid)
call nf95_put_var(ncid, varid_lon, longitude)
call nf95_put_var(ncid, varid_lat, latitude)
call nf95_put_var(ncid, varid_extr_map, extr_map)
call nf95_close(ncid)
print *, 'Created file "extr_map.nc".'
end subroutine write_extr_map
end module write_extr_map_m