! $Id$ module regr_pr_int_m ! Author: Lionel GUEZ implicit none contains SUBROUTINE regr_pr_int(ncid, name, julien, plev, pplay, top_value, v3) ! "regr_pr_int" stands for "regrid pressure interpolation". ! In this procedure: ! -- the root process reads a 2D latitude-pressure field from a ! NetCDF file, at a given day. ! -- the field is packed to the LMDZ horizontal "physics" ! grid and scattered to all threads of all processes; ! -- in all the threads of all the processes, the field is regridded in ! pressure to the LMDZ vertical grid. ! We assume that, in the input file, the field has 3 dimensions: ! latitude, pressure, julian day. ! We assume that latitudes are in ascending order in the input file. ! The target vertical LMDZ grid is the grid of mid-layers. ! Regridding is by linear interpolation. use dimphy, only: klon use netcdf95, only: nf95_inq_varid, nf95_get_var use assert_m, only: assert use regr_lint_m, only: regr_lint use mod_phys_lmdz_mpi_data, only: is_mpi_root use mod_grid_phy_lmdz, only: nbp_lon, nbp_lat, nbp_lev use mod_phys_lmdz_transfert_para, only: scatter2d ! (pack to the LMDZ horizontal "physics" grid and scatter) integer, intent(in):: ncid ! NetCDF ID of the file character(len=*), intent(in):: name ! of the NetCDF variable integer, intent(in):: julien ! jour julien, 1 <= julien <= 360 real, intent(in):: plev(:) ! (pressure level of input data, in Pa, in strictly ascending order) real, intent(in):: pplay(:, :) ! (klon, nbp_lev) ! (pression pour le mileu de chaque couche, en Pa) real, intent(in):: top_value ! (extra value of field at 0 pressure) real, intent(out):: v3(:, :) ! (klon, nbp_lev) ! (regridded field on the partial "physics" grid) ! ("v3(i, k)" is at longitude "xlon(i)", latitude ! "xlat(i)", middle of layer "k".) ! Variables local to the procedure: integer varid, ncerr ! for NetCDF real v1(nbp_lon, nbp_lat, 0:size(plev)) ! (input field at day "julien", on the global "dynamics" horizontal grid) ! (First dimension is for longitude. ! The value is the same for all longitudes. ! "v1(:, j, k >=1)" is at latitude "rlatu(j)" and pressure "plev(k)".) real v2(klon, 0:size(plev)) ! (field scattered to the partial "physics" horizontal grid) ! "v2(i, k >= 1)" is at longitude "xlon(i)", latitude "xlat(i)" ! and pressure "plev(k)".) integer i !-------------------------------------------- CALL assert(shape(v3) == (/klon, nbp_lev/), "regr_pr_int v3") CALL assert(shape(pplay) == (/klon, nbp_lev/), "regr_pr_int pplay") !$omp master if (is_mpi_root) then CALL nf95_inq_varid(ncid, name, varid) ! Get data at the right day from the input file: CALL nf95_get_var(ncid, varid, v1(1, :, 1:), start=(/1, 1, julien/)) ! Latitudes are in ascending order in the input file while ! "rlatu" is in descending order so we need to invert order: v1(1, :, 1:) = v1(1, nbp_lat:1:-1, 1:) ! Complete "v1" with the value at 0 pressure: v1(1, :, 0) = top_value ! Duplicate on all longitudes: v1(2:, :, :) = spread(v1(1, :, :), dim=1, ncopies=nbp_lon-1) end if !$omp end master CALL scatter2d(v1, v2) ! Regrid in pressure at each horizontal position: do i = 1, klon CALL regr_lint(1, v2(i, :), (/0., plev/), pplay(i, nbp_lev:1:-1), & v3(i, nbp_lev:1:-1)) ! (invert order of indices because "pplay" is in descending order) END DO END SUBROUTINE regr_pr_int end module regr_pr_int_m