regr_pr_int_m.f90 Source File


This file depends on

sourcefile~~regr_pr_int_m.f90~~EfferentGraph sourcefile~regr_pr_int_m.f90 regr_pr_int_m.f90 sourcefile~dimphy.f90 dimphy.f90 sourcefile~regr_pr_int_m.f90->sourcefile~dimphy.f90 sourcefile~regr_lint_m.f90 regr_lint_m.f90 sourcefile~regr_pr_int_m.f90->sourcefile~regr_lint_m.f90 sourcefile~mod_phys_lmdz_mpi_data.f90 mod_phys_lmdz_mpi_data.f90 sourcefile~regr_pr_int_m.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_grid_phy_lmdz.f90 mod_grid_phy_lmdz.f90 sourcefile~regr_pr_int_m.f90->sourcefile~mod_grid_phy_lmdz.f90 sourcefile~mod_phys_lmdz_transfert_para.f90 mod_phys_lmdz_transfert_para.f90 sourcefile~regr_pr_int_m.f90->sourcefile~mod_phys_lmdz_transfert_para.f90 sourcefile~assert_m.f90 assert_m.f90 sourcefile~regr_pr_int_m.f90->sourcefile~assert_m.f90 sourcefile~regr_lint_m.f90->sourcefile~assert_m.f90 sourcefile~interpolation.f90 interpolation.f90 sourcefile~regr_lint_m.f90->sourcefile~interpolation.f90 sourcefile~assert_eq_m.f90 assert_eq_m.f90 sourcefile~regr_lint_m.f90->sourcefile~assert_eq_m.f90 sourcefile~lmdz_mpi.f90 lmdz_mpi.F90 sourcefile~mod_phys_lmdz_mpi_data.f90->sourcefile~lmdz_mpi.f90 sourcefile~print_control_mod.f90 print_control_mod.f90 sourcefile~mod_phys_lmdz_mpi_data.f90->sourcefile~print_control_mod.f90 sourcefile~lmdz_cppkeys_wrapper.f90 lmdz_cppkeys_wrapper.F90 sourcefile~mod_phys_lmdz_mpi_data.f90->sourcefile~lmdz_cppkeys_wrapper.f90 sourcefile~mod_phys_lmdz_transfert_para.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_omp_transfert.f90 mod_phys_lmdz_omp_transfert.f90 sourcefile~mod_phys_lmdz_transfert_para.f90->sourcefile~mod_phys_lmdz_omp_transfert.f90 sourcefile~mod_phys_lmdz_mpi_transfert.f90 mod_phys_lmdz_mpi_transfert.f90 sourcefile~mod_phys_lmdz_transfert_para.f90->sourcefile~mod_phys_lmdz_mpi_transfert.f90 sourcefile~mod_phys_lmdz_omp_transfert.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_omp_data.f90 mod_phys_lmdz_omp_data.F90 sourcefile~mod_phys_lmdz_omp_transfert.f90->sourcefile~mod_phys_lmdz_omp_data.f90 sourcefile~mod_phys_lmdz_mpi_transfert.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_mpi_transfert.f90->sourcefile~mod_grid_phy_lmdz.f90 sourcefile~mod_phys_lmdz_mpi_transfert.f90->sourcefile~lmdz_mpi.f90 sourcefile~mod_phys_lmdz_omp_data.f90->sourcefile~dimphy.f90 sourcefile~mod_phys_lmdz_omp_data.f90->sourcefile~mod_phys_lmdz_mpi_data.f90 sourcefile~mod_phys_lmdz_omp_data.f90->sourcefile~print_control_mod.f90

Files dependent on this one

sourcefile~~regr_pr_int_m.f90~~AfferentGraph sourcefile~regr_pr_int_m.f90 regr_pr_int_m.f90 sourcefile~regr_pr_comb_coefoz_m.f90 regr_pr_comb_coefoz_m.f90 sourcefile~regr_pr_comb_coefoz_m.f90->sourcefile~regr_pr_int_m.f90 sourcefile~regr_pr_comb_coefoz_m.f90~2 regr_pr_comb_coefoz_m.f90 sourcefile~regr_pr_comb_coefoz_m.f90~2->sourcefile~regr_pr_int_m.f90 sourcefile~traclmdz_mod.f90 traclmdz_mod.f90 sourcefile~traclmdz_mod.f90->sourcefile~regr_pr_comb_coefoz_m.f90 sourcefile~o3_chem_m.f90 o3_chem_m.f90 sourcefile~traclmdz_mod.f90->sourcefile~o3_chem_m.f90 sourcefile~o3_chem_m.f90->sourcefile~regr_pr_comb_coefoz_m.f90 sourcefile~o3_chem_m.f90~2 o3_chem_m.f90 sourcefile~o3_chem_m.f90~2->sourcefile~regr_pr_comb_coefoz_m.f90 sourcefile~traclmdz_mod.f90~2 traclmdz_mod.f90 sourcefile~traclmdz_mod.f90~2->sourcefile~regr_pr_comb_coefoz_m.f90 sourcefile~traclmdz_mod.f90~2->sourcefile~o3_chem_m.f90 sourcefile~lsc_scav.f90 lsc_scav.f90 sourcefile~lsc_scav.f90->sourcefile~traclmdz_mod.f90 sourcefile~initrrnpb.f90~2 initrrnpb.f90 sourcefile~initrrnpb.f90~2->sourcefile~traclmdz_mod.f90 sourcefile~lsc_scav.f90~2 lsc_scav.f90 sourcefile~lsc_scav.f90~2->sourcefile~traclmdz_mod.f90 sourcefile~lsc_scav_spl.f90 lsc_scav_spl.f90 sourcefile~lsc_scav_spl.f90->sourcefile~traclmdz_mod.f90 sourcefile~cltracrn.f90 cltracrn.f90 sourcefile~cltracrn.f90->sourcefile~traclmdz_mod.f90 sourcefile~initrrnpb.f90 initrrnpb.f90 sourcefile~initrrnpb.f90->sourcefile~traclmdz_mod.f90 sourcefile~lsc_scav_orig.f90 lsc_scav_orig.f90 sourcefile~lsc_scav_orig.f90->sourcefile~traclmdz_mod.f90 sourcefile~radio_decay.f90~2 radio_decay.f90 sourcefile~radio_decay.f90~2->sourcefile~traclmdz_mod.f90 sourcefile~phyredem.f90 phyredem.f90 sourcefile~phyredem.f90->sourcefile~traclmdz_mod.f90 sourcefile~lsc_scav_orig.f90~2 lsc_scav_orig.f90 sourcefile~lsc_scav_orig.f90~2->sourcefile~traclmdz_mod.f90 sourcefile~phyetat0_mod.f90 phyetat0_mod.F90 sourcefile~phyetat0_mod.f90->sourcefile~traclmdz_mod.f90 sourcefile~phytrac_mod.f90 phytrac_mod.f90 sourcefile~phytrac_mod.f90->sourcefile~traclmdz_mod.f90 sourcefile~cltracrn.f90~2 cltracrn.f90 sourcefile~cltracrn.f90~2->sourcefile~traclmdz_mod.f90 sourcefile~lsc_scav_spl.f90~2 lsc_scav_spl.f90 sourcefile~lsc_scav_spl.f90~2->sourcefile~traclmdz_mod.f90 sourcefile~phytrac_mod.f90~2 phytrac_mod.f90 sourcefile~phytrac_mod.f90~2->sourcefile~traclmdz_mod.f90 sourcefile~radio_decay.f90 radio_decay.f90 sourcefile~radio_decay.f90->sourcefile~traclmdz_mod.f90 sourcefile~phys_output_write_mod.f90 phys_output_write_mod.F90 sourcefile~phys_output_write_mod.f90->sourcefile~phytrac_mod.f90 sourcefile~phys_output_write_mod.f90~2 phys_output_write_mod.F90 sourcefile~phys_output_write_mod.f90~2->sourcefile~phytrac_mod.f90 sourcefile~physiqex_mod.f90 physiqex_mod.F90 sourcefile~physiqex_mod.f90->sourcefile~phyetat0_mod.f90 sourcefile~physiq_mod.f90 physiq_mod.F90 sourcefile~physiq_mod.f90->sourcefile~phyetat0_mod.f90 sourcefile~physiq_mod.f90->sourcefile~phytrac_mod.f90 sourcefile~physiq_mod.f90->sourcefile~phys_output_write_mod.f90 sourcefile~physiq_mod.f90->sourcefile~physiqex_mod.f90 sourcefile~diag_slp.f90 diag_slp.f90 sourcefile~physiq_mod.f90->sourcefile~diag_slp.f90 sourcefile~phys_output_mod.f90 phys_output_mod.F90 sourcefile~physiq_mod.f90->sourcefile~phys_output_mod.f90 sourcefile~physiq_mod.f90~2 physiq_mod.F90 sourcefile~physiq_mod.f90~2->sourcefile~phyetat0_mod.f90 sourcefile~physiq_mod.f90~2->sourcefile~phytrac_mod.f90 sourcefile~physiq_mod.f90~2->sourcefile~phys_output_write_mod.f90 sourcefile~physiq_mod.f90~2->sourcefile~physiqex_mod.f90 sourcefile~physiq_mod.f90~2->sourcefile~diag_slp.f90 sourcefile~physiq_mod.f90~2->sourcefile~phys_output_mod.f90 sourcefile~physiqex_mod.f90~2 physiqex_mod.F90 sourcefile~physiqex_mod.f90~2->sourcefile~phyetat0_mod.f90 sourcefile~diag_slp.f90->sourcefile~phys_output_write_mod.f90 sourcefile~phys_output_mod.f90->sourcefile~phys_output_write_mod.f90 sourcefile~old_lmdz1d.f90 old_lmdz1d.f90 sourcefile~old_lmdz1d.f90->sourcefile~physiq_mod.f90 sourcefile~diag_slp.f90~2 diag_slp.f90 sourcefile~diag_slp.f90~2->sourcefile~phys_output_write_mod.f90 sourcefile~callphysiq_mod.f90 callphysiq_mod.f90 sourcefile~callphysiq_mod.f90->sourcefile~physiq_mod.f90 sourcefile~callphysiq_mod.f90~2 callphysiq_mod.f90 sourcefile~callphysiq_mod.f90~2->sourcefile~physiq_mod.f90 sourcefile~scm.f90 scm.f90 sourcefile~scm.f90->sourcefile~physiq_mod.f90 sourcefile~phys_output_mod.f90~2 phys_output_mod.F90 sourcefile~phys_output_mod.f90~2->sourcefile~phys_output_write_mod.f90 sourcefile~recmwf_aero.f90 recmwf_aero.F90 sourcefile~recmwf_aero.f90->sourcefile~phys_output_mod.f90 sourcefile~recmwf_aero.f90~2 recmwf_aero.F90 sourcefile~recmwf_aero.f90~2->sourcefile~phys_output_mod.f90 sourcefile~sw_aeroar4.f90~2 sw_aeroAR4.f90 sourcefile~sw_aeroar4.f90~2->sourcefile~phys_output_mod.f90 sourcefile~calfis.f90 calfis.f90 sourcefile~calfis.f90->sourcefile~callphysiq_mod.f90 sourcefile~sw_aeroar4.f90 sw_aeroAR4.f90 sourcefile~sw_aeroar4.f90->sourcefile~phys_output_mod.f90

Contents

Source Code


Source Code

! $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