source: LMDZ6/trunk/libf/phylmd/regr_pr_int_m.F90 @ 3654

Last change on this file since 3654 was 2788, checked in by dcugnet, 8 years ago

Changes in ce0l about the way ozone forcing files are generated:

1) 3D raw input files "climoz.nc" are now handled.
2) Default behaviour is now to let the gcm interpolate in time online.

This helps to avoid huge forcing files (in particular for 3D fields).
In this case, the output files "climoz_LMDZ.nc" all have 14 records:

  • records 2-13 are obtained with records 1-12 of "climoz.nc".
  • records 1 and 14 are obtained respectively with:
    • record 12 of "climoz_m.nc" if available, of "climoz.nc" otherwise.
    • record 1 of "climoz_p.nc" if available, of "climoz.nc" otherwise.

3) If ok_daily_climoz key is TRUE, the time interpolation (one record

a day) is forced, using the 14 records described below.
This now depends on the calendar (it was on a 360 days basis only).

Changes in the gcm about the way zone forcing files are read/interpolated:

1) 3D horizontally interpolated "climoz_LMDZ.nc" files are now handled.
2) Daily files (already interpolated in time) are still handled, but their

number of records must match the expected number of days, that depends
on the calendar (records step is no longer 1/360 year).

3) 14 records monthly files are now handled (and prefered). This reduces

the I/O to a minimum and the aditional computational cost is low (simple
online linear time interpolation).

4) If adjust_tropopause key is TRUE, the input fields are stretched using

following method:

  • LMDZ dynamical tropopause is detected: Ptrop_lmdz = MAX ( P(Potential Vorticity==2PVU), P(theta==380K) )
  • file chemical tropopause is detected: Ptrop_file = P( tro3 == o3t ), where:

o3t = 91. + 28. * SIN(PI*(month-2)/6) (ppbV)

This formula comes from Thouret & al., ACP 6, 1033-1051, 2006.
The second term of the expression is multiplied by TANH(lat_deg/20.)
to account for latitude dependency.

  • File profile is streched in a +/- 5kms zone around the mean tropopause to ensure resulting tropopause matches the one of LMDZ. See procedure regr_pr_time_av for more details.
  • Property copyright set to
    Name of program: LMDZ
    Creation date: 1984
    Version: LMDZ5
    License: CeCILL version 2
    Holder: Laboratoire de m\'et\'eorologie dynamique, CNRS, UMR 8539
    See the license file in the root directory
File size: 3.7 KB
Line 
1! $Id$
2module regr_pr_int_m
3
4  ! Author: Lionel GUEZ
5
6  implicit none
7
8contains
9
10  subroutine regr_pr_int(ncid, name, julien, plev, pplay, top_value, v3)
11
12    ! "regr_pr_int" stands for "regrid pressure interpolation".
13    ! In this procedure:
14    ! -- the root process reads a 2D latitude-pressure field from a
15    !    NetCDF file, at a given day.
16    ! -- the field is packed to the LMDZ horizontal "physics"
17    !    grid and scattered to all threads of all processes;
18    ! -- in all the threads of all the processes, the field is regridded in
19    !    pressure to the LMDZ vertical grid.
20    ! We assume that, in the input file, the field has 3 dimensions:
21    ! latitude, pressure, julian day.
22    ! We assume that latitudes are in ascending order in the input file.
23    ! The target vertical LMDZ grid is the grid of mid-layers.
24    ! Regridding is by linear interpolation.
25
26    use dimphy, only: klon
27    use netcdf95, only: nf95_inq_varid, handle_err
28    use netcdf, only: nf90_get_var
29    use assert_m, only: assert
30    use regr_lint_m, only: regr_lint
31    use mod_phys_lmdz_mpi_data, only: is_mpi_root
32    use mod_grid_phy_lmdz, only: nbp_lon, nbp_lat, nbp_lev
33    use mod_phys_lmdz_transfert_para, only: scatter2d
34    ! (pack to the LMDZ horizontal "physics" grid and scatter)
35
36    integer, intent(in):: ncid ! NetCDF ID of the file
37    character(len=*), intent(in):: name ! of the NetCDF variable
38    integer, intent(in):: julien ! jour julien, 1 <= julien <= 360
39
40    real, intent(in):: plev(:)
41    ! (pressure level of input data, in Pa, in strictly ascending order)
42
43    real, intent(in):: pplay(:, :) ! (klon, nbp_lev)
44    ! (pression pour le mileu de chaque couche, en Pa)
45
46    real, intent(in):: top_value
47    ! (extra value of field at 0 pressure)
48
49    real, intent(out):: v3(:, :) ! (klon, nbp_lev)
50    ! (regridded field on the partial "physics" grid)
51    ! ("v3(i, k)" is at longitude "xlon(i)", latitude
52    ! "xlat(i)", middle of layer "k".)
53
54    ! Variables local to the procedure:
55
56    integer varid, ncerr ! for NetCDF
57
58    real  v1(nbp_lon, nbp_lat, 0:size(plev))
59    ! (input field at day "julien", on the global "dynamics" horizontal grid)
60    ! (First dimension is for longitude.
61    ! The value is the same for all longitudes.
62    ! "v1(:, j, k >=1)" is at latitude "rlatu(j)" and pressure "plev(k)".)
63
64    real v2(klon, 0:size(plev))
65    ! (field scattered to the partial "physics" horizontal grid)
66    ! "v2(i, k >= 1)" is at longitude "xlon(i)", latitude "xlat(i)"
67    ! and pressure "plev(k)".)
68
69    integer i
70
71    !--------------------------------------------
72
73    call assert(shape(v3) == (/klon, nbp_lev/), "regr_pr_int v3")
74    call assert(shape(pplay) == (/klon, nbp_lev/), "regr_pr_int pplay")
75
76    !$omp master
77    if (is_mpi_root) then
78       call nf95_inq_varid(ncid, name, varid)
79
80       ! Get data at the right day from the input file:
81       ncerr = nf90_get_var(ncid, varid, v1(1, :, 1:), start=(/1, 1, julien/))
82       call handle_err("regr_pr_int nf90_get_var " // name, ncerr, ncid)
83       ! Latitudes are in ascending order in the input file while
84       ! "rlatu" is in descending order so we need to invert order:
85       v1(1, :, 1:) = v1(1, nbp_lat:1:-1, 1:)
86
87       ! Complete "v1" with the value at 0 pressure:
88       v1(1, :, 0) = top_value
89
90       ! Duplicate on all longitudes:
91       v1(2:, :, :) = spread(v1(1, :, :), dim=1, ncopies=nbp_lon-1)
92    end if
93    !$omp end master
94
95    call scatter2d(v1, v2)
96
97    ! Regrid in pressure at each horizontal position:
98    do i = 1, klon
99       call regr_lint(1, v2(i, :), (/0., plev/), pplay(i, nbp_lev:1:-1), &
100                         v3(i, nbp_lev:1:-1))
101       ! (invert order of indices because "pplay" is in descending order)
102    end do
103
104  end subroutine regr_pr_int
105
106end module regr_pr_int_m
Note: See TracBrowser for help on using the repository browser.