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, handle_err |
---|
28 | use netcdf, only: nf90_nowrite, nf90_get_var |
---|
29 | use assert_m, only: assert |
---|
30 | use regr1_step_av_m, only: regr1_step_av |
---|
31 | use press_coefoz_m, only: press_in_edg |
---|
32 | use control_mod, only: dayref |
---|
33 | |
---|
34 | REAL, intent(in):: p3d(:, :, :) ! pressure at layer interfaces, in Pa |
---|
35 | ! ("p3d(i, j, l)" is at longitude "rlonv(i)", latitude "rlatu(j)", |
---|
36 | ! for interface "l") |
---|
37 | |
---|
38 | real, intent(out):: o3_mob_regr(:, :, :) ! (iim + 1, jjm + 1, llm) |
---|
39 | ! (ozone mole fraction from Mobidic adapted to the LMDZ grid) |
---|
40 | ! ("o3_mob_regr(i, j, l)" is at longitude "rlonv(i)", latitude |
---|
41 | ! "rlatu(j)" and pressure level "pls(i, j, l)") |
---|
42 | |
---|
43 | ! Variables local to the procedure: |
---|
44 | |
---|
45 | include "dimensions.h" |
---|
46 | |
---|
47 | integer ncid, varid, ncerr ! for NetCDF |
---|
48 | integer i, j |
---|
49 | |
---|
50 | real r_mob(jjm + 1, size(press_in_edg) - 1) |
---|
51 | ! (ozone mole fraction from Mobidic at day "dayref") |
---|
52 | ! (r_mob(j, k) is at latitude "rlatu(j)", in pressure interval |
---|
53 | ! "[press_in_edg(k), press_in_edg(k+1)]".) |
---|
54 | |
---|
55 | !------------------------------------------------------------ |
---|
56 | |
---|
57 | print *, "Call sequence information: regr_pr_o3" |
---|
58 | call assert(shape(o3_mob_regr) == (/iim + 1, jjm + 1, llm/), & |
---|
59 | "regr_pr_o3 o3_mob_regr") |
---|
60 | call assert(shape(p3d) == (/iim + 1, jjm + 1, llm + 1/), "regr_pr_o3 p3d") |
---|
61 | |
---|
62 | call nf95_open("coefoz_LMDZ.nc", nf90_nowrite, ncid) |
---|
63 | |
---|
64 | call nf95_inq_varid(ncid, "r_Mob", varid) |
---|
65 | ! Get data at the right day from the input file: |
---|
66 | ncerr = nf90_get_var(ncid, varid, r_mob, start=(/1, 1, dayref/)) |
---|
67 | call handle_err("nf90_get_var r_Mob", ncerr) |
---|
68 | ! Latitudes are in ascending order in the input file while |
---|
69 | ! "rlatu" is in descending order so we need to invert order: |
---|
70 | r_mob = r_mob(jjm+1:1:-1, :) |
---|
71 | |
---|
72 | call nf95_close(ncid) |
---|
73 | |
---|
74 | ! Regrid in pressure by averaging a step function of pressure: |
---|
75 | |
---|
76 | ! Poles: |
---|
77 | do j = 1, jjm + 1, jjm |
---|
78 | o3_mob_regr(1, j, llm:1:-1) & |
---|
79 | = regr1_step_av(r_mob(j, :), press_in_edg, p3d(1, j, llm+1:1:-1)) |
---|
80 | ! (invert order of indices because "p3d" is in descending order) |
---|
81 | end do |
---|
82 | |
---|
83 | ! Other latitudes: |
---|
84 | do j = 2, jjm |
---|
85 | do i = 1, iim |
---|
86 | o3_mob_regr(i, j, llm:1:-1) & |
---|
87 | = regr1_step_av(r_mob(j, :), press_in_edg, & |
---|
88 | p3d(i, j, llm+1:1:-1)) |
---|
89 | ! (invert order of indices because "p3d" is in descending order) |
---|
90 | end do |
---|
91 | end do |
---|
92 | |
---|
93 | ! Duplicate pole values on all longitudes: |
---|
94 | o3_mob_regr(2:, 1, :) = spread(o3_mob_regr(1, 1, :), dim=1, ncopies=iim) |
---|
95 | o3_mob_regr(2:, jjm + 1, :) & |
---|
96 | = spread(o3_mob_regr(1, jjm + 1, :), dim=1, ncopies=iim) |
---|
97 | |
---|
98 | ! Duplicate first longitude to last longitude: |
---|
99 | o3_mob_regr(iim + 1, 2:jjm, :) = o3_mob_regr(1, 2:jjm, :) |
---|
100 | |
---|
101 | end subroutine regr_pr_o3 |
---|
102 | |
---|
103 | end module regr_pr_o3_m |
---|