Ignore:
Timestamp:
Feb 9, 2016, 3:45:31 PM (8 years ago)
Author:
lguez
Message:

For read_climoz = 1 or 2, replaced first order conservative regridding
of ozone by second order conservative regridding, with Van Leer
slope-limiting. The replacement is done for both latitude and pressure
regridding. The replacement is beneficial if the resolution of the
input data is coarser than the resolution of LMDZ. If the resolution
of the input data is finer, then the replacement is neutral, it does
not change much.

File:
1 edited

Legend:

Unmodified
Added
Removed
  • LMDZ5/trunk/libf/phylmd/regr_lat_time_climoz_m.F90

    r2346 r2440  
    6767
    6868    use mod_grid_phy_lmdz, ONLY : nbp_lat
    69     use regr1_step_av_m, only: regr1_step_av
     69    use regr1_conserv_m, only: regr1_conserv
    7070    use regr3_lint_m, only: regr3_lint
    7171    use netcdf95, only: handle_err, nf95_close, nf95_get_att, nf95_gw_var, &
     
    7676    use regular_lonlat_mod, only : boundslat_reg, south
    7777    use nrtype, only: pi
     78    use slopes_m, only: slopes
    7879
    7980    integer, intent(in):: read_climoz ! read ozone climatology
     
    9293    ! (of input data, converted to rad, sorted in strictly ascending order)
    9394
    94     real, allocatable:: lat_in_edg(:)
    95     ! (edges of latitude intervals for input data, in rad, in strictly
     95    real, allocatable:: sin_lat_in_edg(:)
     96    ! (sine of edges of latitude intervals for input data, in rad, in strictly
    9697    ! ascending order)
    9798
     
    115116
    116117    real, allocatable:: o3_regr_lat(:, :, :, :)
    117     ! (jjm + 1, n_plev, 0:13, read_climoz)
     118    ! (nbp_lat, n_plev, 0:13, read_climoz)
    118119    ! mean of "o3_in" over a latitude interval of LMDZ
    119120    ! First dimension is latitude interval.
    120121    ! The latitude interval for "o3_regr_lat(j,:, :, :)" contains "rlatu(j)".
    121     ! If "j" is between 2 and "jjm" then the interval is:
     122    ! If "j" is between 2 and "nbp_lat - 1" then the interval is:
    122123    ! [rlatv(j), rlatv(j-1)]
    123     ! If "j" is 1 or "jjm + 1" then the interval is:
     124    ! If "j" is 1 or "nbp_lat" then the interval is:
    124125    ! [rlatv(1), pi / 2]
    125126    ! or:
    126     ! [- pi / 2, rlatv(jjm)]
     127    ! [- pi / 2, rlatv(nbp_lat - 1)]
    127128    ! respectively.
    128129    ! "o3_regr_lat(:, k, :, :)" is for pressure level "plev(k)".
     
    132133
    133134    real, allocatable:: o3_out(:, :, :, :)
    134     ! (jjm + 1, n_plev, 360, read_climoz)
     135    ! (nbp_lat, n_plev, 360, read_climoz)
    135136    ! regridded ozone climatology
    136137    ! "o3_out(j, k, l, :)" is at latitude "rlatu(j)", pressure
     
    175176    latitude = latitude / 180. * pi
    176177    n_lat = size(latitude)
    177     ! We need to supply the latitudes to "regr1_step_av" in
     178    ! We need to supply the latitudes to "regr1_conserv" in
    178179    ! ascending order, so invert order if necessary:
    179180    desc_lat = latitude(1) > latitude(n_lat)
     
    181182
    182183    ! Compute edges of latitude intervals:
    183     allocate(lat_in_edg(n_lat + 1))
    184     lat_in_edg(1) = - pi / 2
    185     forall (j = 2:n_lat) lat_in_edg(j) = (latitude(j - 1) + latitude(j)) / 2
    186     lat_in_edg(n_lat + 1) = pi / 2
     184    allocate(sin_lat_in_edg(n_lat + 1))
     185    sin_lat_in_edg(1) = - 1.
     186    forall (j = 2:n_lat) sin_lat_in_edg(j) = sin((latitude(j - 1) &
     187         + latitude(j)) / 2.)
     188    sin_lat_in_edg(n_lat + 1) = 1.
    187189    deallocate(latitude) ! pointer
    188190
     
    292294       print *, &
    293295            "Found 12 months in ozone climatologies, assuming periodicity..."
    294        o3_regr_lat(nbp_lat:1:-1, :, 1:12, :) = regr1_step_av(o3_in, &
    295             xs=sin(lat_in_edg), xt=sin((/- pi / 2, boundslat_reg(nbp_lat-1:1:-1,south), pi / 2/)))
     296       call regr1_conserv(o3_in, xs = sin_lat_in_edg, &
     297            xt = (/- 1., sin(boundslat_reg(nbp_lat - 1:1:- 1, south)), 1./), &
     298            vt = o3_regr_lat(nbp_lat:1:- 1, :, 1:12, :), &
     299            slope = slopes(o3_in, sin_lat_in_edg))
    296300       ! (invert order of indices in "o3_regr_lat" because "rlatu" is
    297301       ! in descending order)
     
    303307    else
    304308       print *, "Using 14 months in ozone climatologies..."
    305        o3_regr_lat(nbp_lat:1:-1, :, :, :) = regr1_step_av(o3_in, &
    306             xs=sin(lat_in_edg), xt=sin((/- pi / 2, boundslat_reg(nbp_lat-1:1:-1,south), pi / 2/)))
     309       call regr1_conserv(o3_in, xs = sin_lat_in_edg, &
     310            xt = (/- 1., sin(boundslat_reg(nbp_lat - 1:1:- 1, south)), 1./), &
     311            vt = o3_regr_lat(nbp_lat:1:- 1, :, :, :), &
     312            slope = slopes(o3_in, sin_lat_in_edg))
    307313       ! (invert order of indices in "o3_regr_lat" because "rlatu" is
    308314       ! in descending order)
Note: See TracChangeset for help on using the changeset viewer.