/[lmdze]/trunk/libf/phylmd/Mobidic/regr_pr_coefoz.f90
ViewVC logotype

Annotation of /trunk/libf/phylmd/Mobidic/regr_pr_coefoz.f90

Parent Directory Parent Directory | Revision Log Revision Log


Revision 17 - (hide annotations)
Tue Aug 5 13:31:32 2008 UTC (15 years, 9 months ago) by guez
File size: 4846 byte(s)
Created rule for "compare_sampl_*" files in
"Documentation/Manuel_LMDZE.texfol/Graphiques/GNUmakefile".

Extracted "qcheck", "radiornpb", "minmaxqfi" into separate files.

Read pressure coordinate of ozone coefficients once per run instead of
every day.

Added some "intent" attributes.

Added argument "nq" to "ini_histday". Replaced calls to "gr_fi_ecrit"
by calls to "gr_phy_write_2d". "Sigma_O3_Royer" is written to
"histday.nc" only if "nq >= 4". Moved "ini_histrac" to module
"ini_hist".

Compute "zmasse" in "physiq", pass it to "phytrac".

Removed computations of "pftsol*" and "ppsrf*" in "phytrac".

Do not use variable "rg" from module "YOMCST" in "TLIFT".

1 guez 7 module regr_pr_coefoz
2 guez 3
3     implicit none
4    
5     contains
6    
7 guez 17 subroutine regr_pr_av_coefoz(ncid, name, julien, v3)
8 guez 3
9 guez 7 ! "regr_pr_av_coefoz" stands for "regrid pressure averaging
10 guez 10 ! coefficient ozone".
11 guez 7 ! This procedure reads a single Mobidic ozone coefficient from
12 guez 12 ! "coefoz_LMDZ.nc", at the current day, regrids this parameter in
13 guez 7 ! pressure to the LMDZ vertical grid and packs it to the LMDZ
14     ! horizontal "physics" grid.
15 guez 10 ! Regridding in pressure is done by averaging a step function.
16 guez 3
17     use dimens_m, only: iim, jjm, llm
18     use dimphy, only: klon
19 guez 7 use netcdf95, only: nf95_inq_varid, handle_err
20     use netcdf, only: nf90_get_var
21 guez 3 use grid_change, only: dyn_phy
22 guez 7 use regr_pr, only: regr_pr_av
23 guez 13 use numer_rec, only: assert
24 guez 17 use press_coefoz_m, only: press_in_edg
25 guez 3
26     integer, intent(in):: ncid ! NetCDF ID of the file
27     character(len=*), intent(in):: name ! of the NetCDF variable
28 guez 7 integer, intent(in):: julien ! jour julien, 1 <= julien <= 360
29 guez 3
30 guez 7 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 guez 17 subroutine regr_pr_int_coefoz(ncid, name, julien, top_value, v3)
72 guez 7
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 guez 13 use numer_rec, only: assert
86 guez 17 use press_coefoz_m, only: plev
87 guez 7
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 guez 3 ! (ozone parameter from Mobidic on the "physics" grid)
97 guez 7 ! ("v3(i, k)" is at longitude "xlon(i)", latitude
98 guez 3 ! "xlat(i)", middle of layer "k".)
99    
100     ! Variables local to the procedure:
101     integer varid, ncerr
102 guez 7 integer k
103 guez 3
104 guez 7 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 guez 3 ! (ozone parameter from Mobidic on the "dynamics" grid)
110 guez 7 ! "v2(i, j, k)" is at longitude "rlonv(i)", latitude
111 guez 3 ! "rlatu(j)", middle of layer "k".)
112    
113     !--------------------------------------------
114    
115 guez 7 call assert(shape(v3) == (/klon, llm/), "regr_pr_int_coefoz")
116    
117 guez 3 call nf95_inq_varid(ncid, name, varid)
118    
119 guez 7 ! 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 guez 3 ! Latitudes are in increasing order in the input file while
123 guez 7 ! "rlatu" is in decreasing order so we need to invert order:
124     v1(:, 1:) = v1(jjm+1:1:-1, 1:)
125 guez 3
126 guez 7 ! Complete "v1" with the value at 0 pressure:
127     v1(:, 0) = top_value
128 guez 3
129 guez 7 ! 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

  ViewVC Help
Powered by ViewVC 1.1.21