1 |
guez |
8 |
module regr_pr_o3_m |
2 |
|
|
|
3 |
|
|
implicit none |
4 |
|
|
|
5 |
|
|
contains |
6 |
|
|
|
7 |
|
|
subroutine regr_pr_o3(o3_mob_regr) |
8 |
|
|
|
9 |
|
|
! "regr_pr_o3" stands for "regrid pressure ozone". |
10 |
|
|
! This procedure reads Mobidic ozone mole fraction from |
11 |
|
|
! "coefoz_LMDZ.nc" at the initial day and regrids it in pressure. |
12 |
|
|
! Regridding is by averaging, assuming a step function. |
13 |
|
|
! We assume that, in the input file, the pressure levels are in hPa |
14 |
|
|
! and strictly increasing. |
15 |
|
|
|
16 |
|
|
use conf_gcm_m, only: dayref |
17 |
|
|
use dimens_m, only: iim, jjm, llm |
18 |
|
|
use netcdf95, only: nf95_open, nf95_close, nf95_inq_varid, handle_err, & |
19 |
|
|
nf95_get_coord |
20 |
|
|
use netcdf, only: nf90_nowrite, nf90_get_var |
21 |
|
|
use regr_pr, only: regr_pr_av |
22 |
guez |
13 |
use numer_rec, only: assert |
23 |
guez |
8 |
|
24 |
|
|
real, intent(out):: o3_mob_regr(:, :, :) ! (iim + 1, jjm + 1, llm) |
25 |
|
|
! (ozone mole fraction from Mobidic adapted to the LMDZ grid) |
26 |
|
|
! ("o3_mob_regr(i, j, l)" is at longitude "rlonv(i)", latitude |
27 |
|
|
! "rlatu(j)" and pressure level "pls(i, j, l)") |
28 |
|
|
|
29 |
|
|
! Variables local to the procedure: |
30 |
|
|
|
31 |
|
|
real, pointer:: plev(:) |
32 |
|
|
! (pressure levels of Mobidic data, in Pa, in strictly increasing order) |
33 |
|
|
|
34 |
|
|
real, allocatable:: press_in_edg(:) |
35 |
|
|
! (edges of pressure intervals for Mobidic data, in Pa, in strictly |
36 |
|
|
! increasing order) |
37 |
|
|
|
38 |
|
|
integer ncid, varid, ncerr ! for NetCDF |
39 |
|
|
integer n_plev ! number of pressure levels in Mobidic data |
40 |
|
|
integer j |
41 |
|
|
|
42 |
|
|
real, allocatable:: r_mob(:, :)! (jjm + 1, n_plev) |
43 |
|
|
! (ozone mole fraction from Mobidic at day "dayref") |
44 |
|
|
! (r_mob(j, k) is at latitude "rlatu(j)" and pressure level "plev(k)".) |
45 |
|
|
|
46 |
|
|
!------------------------------------------------------------ |
47 |
|
|
|
48 |
guez |
10 |
print *, "Call sequence information: regr_pr_o3" |
49 |
guez |
8 |
call assert(shape(o3_mob_regr) == (/iim + 1, jjm + 1, llm/), "regr_pr_o3") |
50 |
|
|
|
51 |
|
|
call nf95_open("coefoz_LMDZ.nc", nf90_nowrite, ncid) |
52 |
|
|
|
53 |
|
|
call nf95_get_coord(ncid, "plev", plev) |
54 |
|
|
! Convert from hPa to Pa because "regr_pr_av" requires so: |
55 |
|
|
plev = plev * 100. |
56 |
|
|
n_plev = size(plev) |
57 |
|
|
|
58 |
|
|
! Compute edges of pressure intervals: |
59 |
|
|
allocate(press_in_edg(n_plev + 1)) |
60 |
|
|
press_in_edg(1) = 0. |
61 |
|
|
! We choose edges halfway in logarithm: |
62 |
|
|
forall (j = 2:n_plev) press_in_edg(j) = sqrt(plev(j - 1) * plev(j)) |
63 |
|
|
press_in_edg(n_plev + 1) = huge(0.) |
64 |
|
|
! (infinity, but any value guaranteed to be greater than the |
65 |
|
|
! surface pressure would do) |
66 |
|
|
|
67 |
|
|
deallocate(plev) ! pointer |
68 |
|
|
|
69 |
|
|
call nf95_inq_varid(ncid, "r_Mob", varid) |
70 |
|
|
allocate(r_mob(jjm + 1, n_plev)) |
71 |
|
|
|
72 |
|
|
! Get data at the right day from the input file: |
73 |
|
|
ncerr = nf90_get_var(ncid, varid, r_mob, start=(/1, 1, dayref/)) |
74 |
|
|
call handle_err("nf90_get_var r_Mob", ncerr) |
75 |
|
|
! Latitudes are in increasing order in the input file while |
76 |
|
|
! "rlatu" is in decreasing order so we need to invert order: |
77 |
|
|
r_mob = r_mob(jjm+1:1:-1, :) |
78 |
|
|
|
79 |
|
|
call nf95_close(ncid) |
80 |
|
|
|
81 |
|
|
o3_mob_regr = regr_pr_av(r_mob, press_in_edg) |
82 |
|
|
|
83 |
|
|
end subroutine regr_pr_o3 |
84 |
|
|
|
85 |
|
|
end module regr_pr_o3_m |