1 |
|
|
! $Id$ |
2 |
|
|
module regr_pr_o3_m |
3 |
|
|
|
4 |
|
|
implicit none |
5 |
|
|
|
6 |
|
|
contains |
7 |
|
|
|
8 |
|
|
subroutine regr_pr_o3(p3d, o3_mob_regr) |
9 |
|
|
|
10 |
|
|
! "regr_pr_o3" stands for "regrid pressure ozone". |
11 |
|
|
! This procedure reads Mobidic ozone mole fraction from |
12 |
|
|
! "coefoz_LMDZ.nc" at the initial day of the run and regrids it in |
13 |
|
|
! pressure. |
14 |
|
|
! Ozone mole fraction from "coefoz_LMDZ.nc" at the initial day is |
15 |
|
|
! a 2D latitude -- pressure variable. |
16 |
|
|
! The target horizontal LMDZ grid is the "scalar" grid: "rlonv", "rlatu". |
17 |
|
|
! The target vertical LMDZ grid is the grid of layer boundaries. |
18 |
|
|
! We assume that the input variable is already on the LMDZ "rlatu" |
19 |
|
|
! latitude grid. |
20 |
|
|
! The input variable does not depend on longitude, but the |
21 |
|
|
! pressure at LMDZ layers does. |
22 |
|
|
! Therefore, the values on the LMDZ grid do depend on longitude. |
23 |
|
|
! Regridding is by averaging, assuming a step function. |
24 |
|
|
! We assume that, in the input file, the pressure levels are in |
25 |
|
|
! hPa and strictly increasing. |
26 |
|
|
|
27 |
|
|
use netcdf95, only: nf95_open, nf95_close, nf95_inq_varid, nf95_get_var |
28 |
|
|
use netcdf, only: nf90_nowrite |
29 |
|
|
use assert_m, only: assert |
30 |
|
|
use regr_conserv_m, only: regr_conserv |
31 |
|
|
use press_coefoz_m, only: press_in_edg |
32 |
|
|
use time_phylmdz_mod, only: day_ref |
33 |
|
|
use mod_grid_phy_lmdz, only: nbp_lon, nbp_lat, nbp_lev |
34 |
|
|
|
35 |
|
|
REAL, intent(in):: p3d(:, :, :) ! pressure at layer interfaces, in Pa |
36 |
|
|
! ("p3d(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)", |
37 |
|
|
! for interface "l") |
38 |
|
|
|
39 |
|
|
real, intent(out):: o3_mob_regr(:, :, :) ! (iim + 1, jjm + 1, llm) |
40 |
|
|
! (ozone mole fraction from Mobidic adapted to the LMDZ grid) |
41 |
|
|
! ("o3_mob_regr(i, j, l)" is at longitude "rlonv(i)", latitude |
42 |
|
|
! "rlatu(j)" and pressure level "pls(i, j, l)") |
43 |
|
|
|
44 |
|
|
! Variables local to the procedure: |
45 |
|
|
|
46 |
|
|
integer ncid, varid, ncerr ! for NetCDF |
47 |
|
|
integer i, j |
48 |
|
|
|
49 |
|
|
real r_mob(nbp_lat, size(press_in_edg) - 1) |
50 |
|
|
! (ozone mole fraction from Mobidic at day "day_ref") |
51 |
|
|
! (r_mob(j, k) is at latitude "rlatu(j)", in pressure interval |
52 |
|
|
! "[press_in_edg(k), press_in_edg(k+1)]".) |
53 |
|
|
|
54 |
|
|
!------------------------------------------------------------ |
55 |
|
|
|
56 |
|
|
print *, "Call sequence information: regr_pr_o3" |
57 |
|
|
call assert(shape(o3_mob_regr) == (/nbp_lon + 1, nbp_lat, nbp_lev/), & |
58 |
|
|
"regr_pr_o3 o3_mob_regr") |
59 |
|
|
call assert(shape(p3d) == (/nbp_lon + 1, nbp_lat, nbp_lev + 1/), "regr_pr_o3 p3d") |
60 |
|
|
|
61 |
|
|
call nf95_open("coefoz_LMDZ.nc", nf90_nowrite, ncid) |
62 |
|
|
|
63 |
|
|
call nf95_inq_varid(ncid, "r_Mob", varid) |
64 |
|
|
! Get data at the right day from the input file: |
65 |
|
|
call nf95_get_var(ncid, varid, r_mob, start=(/1, 1, day_ref/)) |
66 |
|
|
! Latitudes are in ascending order in the input file while |
67 |
|
|
! "rlatu" is in descending order so we need to invert order: |
68 |
|
|
r_mob = r_mob(nbp_lat:1:-1, :) |
69 |
|
|
|
70 |
|
|
call nf95_close(ncid) |
71 |
|
|
|
72 |
|
|
! Regrid in pressure by averaging a step function of pressure: |
73 |
|
|
|
74 |
|
|
! Poles: |
75 |
|
|
do j = 1, nbp_lat, nbp_lat-1 |
76 |
|
|
call regr_conserv(1, r_mob(j, :), press_in_edg, & |
77 |
|
|
p3d(1, j, nbp_lev + 1:1:-1), o3_mob_regr(1, j, nbp_lev:1:-1)) |
78 |
|
|
! (invert order of indices because "p3d" is in descending order) |
79 |
|
|
end do |
80 |
|
|
|
81 |
|
|
! Other latitudes: |
82 |
|
|
do j = 2, nbp_lat-1 |
83 |
|
|
do i = 1, nbp_lon |
84 |
|
|
call regr_conserv(1, r_mob(j, :), press_in_edg, & |
85 |
|
|
p3d(i, j, nbp_lev + 1:1:-1), o3_mob_regr(i, j, nbp_lev:1:-1)) |
86 |
|
|
! (invert order of indices because "p3d" is in descending order) |
87 |
|
|
end do |
88 |
|
|
end do |
89 |
|
|
|
90 |
|
|
! Duplicate pole values on all longitudes: |
91 |
|
|
o3_mob_regr(2:, 1, :) = spread(o3_mob_regr(1, 1, :), dim=1, ncopies=nbp_lon) |
92 |
|
|
o3_mob_regr(2:, nbp_lat, :) & |
93 |
|
|
= spread(o3_mob_regr(1, nbp_lat, :), dim=1, ncopies=nbp_lon) |
94 |
|
|
|
95 |
|
|
! Duplicate first longitude to last longitude: |
96 |
|
|
o3_mob_regr(nbp_lon + 1, 2:nbp_lat-1, :) = o3_mob_regr(1, 2:nbp_lat-1, :) |
97 |
|
|
|
98 |
|
|
end subroutine regr_pr_o3 |
99 |
|
|
|
100 |
|
|
end module regr_pr_o3_m |