1 |
module regr_pr_coefoz |
2 |
|
3 |
implicit none |
4 |
|
5 |
contains |
6 |
|
7 |
subroutine regr_pr_av_coefoz(ncid, name, julien, v3) |
8 |
|
9 |
! "regr_pr_av_coefoz" stands for "regrid pressure averaging |
10 |
! coefficient ozone". |
11 |
! This procedure reads a single Mobidic ozone coefficient from |
12 |
! "coefoz_LMDZ.nc", at the current day, regrids this parameter in |
13 |
! pressure to the LMDZ vertical grid and packs it to the LMDZ |
14 |
! horizontal "physics" grid. |
15 |
! Regridding in pressure is done by averaging a step function. |
16 |
|
17 |
use dimens_m, only: iim, jjm, llm |
18 |
use dimphy, only: klon |
19 |
use netcdf95, only: nf95_inq_varid, handle_err |
20 |
use netcdf, only: nf90_get_var |
21 |
use grid_change, only: dyn_phy |
22 |
use regr_pr, only: regr_pr_av |
23 |
use numer_rec, only: assert |
24 |
use press_coefoz_m, only: press_in_edg |
25 |
|
26 |
integer, intent(in):: ncid ! NetCDF ID of the file |
27 |
character(len=*), intent(in):: name ! of the NetCDF variable |
28 |
integer, intent(in):: julien ! jour julien, 1 <= julien <= 360 |
29 |
|
30 |
real, intent(out):: v3(:, :) ! (klon, llm) |
31 |
! (ozone coefficient from Mobidic on the "physics" grid) |
32 |
! ("v3(i, k)" is at longitude "xlon(i)", latitude |
33 |
! "xlat(i)", middle of layer "k".) |
34 |
|
35 |
! Variables local to the procedure: |
36 |
integer varid, ncerr |
37 |
integer k |
38 |
|
39 |
real v1(jjm + 1, size(press_in_edg) - 1) |
40 |
! (ozone coefficient from "coefoz_LMDZ.nc" at day "julien") |
41 |
! ("v1(j, k)" is at latitude "rlatu(j)" and for |
42 |
! pressure interval "[press_in_edg(k), press_in_edg(k+1)]".) |
43 |
|
44 |
real v2(iim + 1, jjm + 1, llm) |
45 |
! (ozone parameter from Mobidic on the "dynamics" grid) |
46 |
! "v2(i, j, k)" is at longitude "rlonv(i)", latitude |
47 |
! "rlatu(j)", middle of layer "k".) |
48 |
|
49 |
!-------------------------------------------- |
50 |
|
51 |
call assert(shape(v3) == (/klon, llm/), "regr_pr_av_coefoz") |
52 |
|
53 |
call nf95_inq_varid(ncid, name, varid) |
54 |
|
55 |
! Get data at the right day from the input file: |
56 |
ncerr = nf90_get_var(ncid, varid, v1, start=(/1, 1, julien/)) |
57 |
call handle_err("regr_pr_av_coefoz nf90_get_var " // name, ncerr, ncid) |
58 |
! Latitudes are in increasing order in the input file while |
59 |
! "rlatu" is in decreasing order so we need to invert order: |
60 |
v1 = v1(jjm+1:1:-1, :) |
61 |
|
62 |
! Regrid in pressure at each horizontal position: |
63 |
v2 = regr_pr_av(v1, press_in_edg) |
64 |
|
65 |
forall (k = 1:llm) v3(:, k) = pack(v2(:, :, k), dyn_phy) |
66 |
|
67 |
end subroutine regr_pr_av_coefoz |
68 |
|
69 |
!*************************************************************** |
70 |
|
71 |
subroutine regr_pr_int_coefoz(ncid, name, julien, top_value, v3) |
72 |
|
73 |
! This procedure reads a single Mobidic ozone coefficient from |
74 |
!"coefoz_LMDZ.nc", at the current day, regrids this parameter in |
75 |
! pressure to the LMDZ vertical grid and packs it to the LMDZ |
76 |
! horizontal "physics" grid. |
77 |
! Regridding is by linear interpolation. |
78 |
|
79 |
use dimens_m, only: iim, jjm, llm |
80 |
use dimphy, only: klon |
81 |
use netcdf95, only: nf95_inq_varid, handle_err |
82 |
use netcdf, only: nf90_get_var |
83 |
use grid_change, only: dyn_phy |
84 |
use regr_pr, only: regr_pr_int |
85 |
use numer_rec, only: assert |
86 |
use press_coefoz_m, only: plev |
87 |
|
88 |
integer, intent(in):: ncid ! NetCDF ID of the file |
89 |
character(len=*), intent(in):: name ! of the NetCDF variable |
90 |
integer, intent(in):: julien ! jour julien, 1 <= julien <= 360 |
91 |
|
92 |
real, intent(in):: top_value |
93 |
! (extra value of ozone coefficient at 0 pressure) |
94 |
|
95 |
real, intent(out):: v3(:, :) ! (klon, llm) |
96 |
! (ozone parameter from Mobidic on the "physics" grid) |
97 |
! ("v3(i, k)" is at longitude "xlon(i)", latitude |
98 |
! "xlat(i)", middle of layer "k".) |
99 |
|
100 |
! Variables local to the procedure: |
101 |
integer varid, ncerr |
102 |
integer k |
103 |
|
104 |
real v1(jjm + 1, 0:size(plev)) |
105 |
! (ozone coefficient from "coefoz_LMDZ.nc" at day "julien") |
106 |
! ("v1(j, k >=1)" is at latitude "rlatu(j)" and pressure "plev(k)".) |
107 |
|
108 |
real v2(iim + 1, jjm + 1, llm) |
109 |
! (ozone parameter from Mobidic on the "dynamics" grid) |
110 |
! "v2(i, j, k)" is at longitude "rlonv(i)", latitude |
111 |
! "rlatu(j)", middle of layer "k".) |
112 |
|
113 |
!-------------------------------------------- |
114 |
|
115 |
call assert(shape(v3) == (/klon, llm/), "regr_pr_int_coefoz") |
116 |
|
117 |
call nf95_inq_varid(ncid, name, varid) |
118 |
|
119 |
! Get data at the right day from the input file: |
120 |
ncerr = nf90_get_var(ncid, varid, v1(:, 1:), start=(/1, 1, julien/)) |
121 |
call handle_err("regr_pr_int_coefoz nf90_get_var " // name, ncerr, ncid) |
122 |
! Latitudes are in increasing order in the input file while |
123 |
! "rlatu" is in decreasing order so we need to invert order: |
124 |
v1(:, 1:) = v1(jjm+1:1:-1, 1:) |
125 |
|
126 |
! Complete "v1" with the value at 0 pressure: |
127 |
v1(:, 0) = top_value |
128 |
|
129 |
! Regrid in pressure at each horizontal position: |
130 |
v2 = regr_pr_int(v1, (/0., plev/)) |
131 |
|
132 |
forall (k = 1:llm) v3(:, k) = pack(v2(:, :, k), dyn_phy) |
133 |
|
134 |
end subroutine regr_pr_int_coefoz |
135 |
|
136 |
end module regr_pr_coefoz |