/[lmdze]/trunk/phylmd/Mobidic/press_coefoz.f
ViewVC logotype

Annotation of /trunk/phylmd/Mobidic/press_coefoz.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 17 - (hide annotations)
Tue Aug 5 13:31:32 2008 UTC (15 years, 10 months ago) by guez
Original Path: trunk/libf/phylmd/Mobidic/press_coefoz.f90
File size: 1568 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 17 module press_coefoz_m
2    
3     implicit none
4    
5     real, pointer, save:: plev(:)
6     ! (pressure level of Mobidic input data, converted to Pa, in
7     ! ascending order)
8    
9     real, allocatable, save:: press_in_edg(:)
10     ! (edges of pressure intervals for Mobidic input data, in Pa, in
11     ! ascending order)
12    
13     contains
14    
15     subroutine press_coefoz
16    
17     ! This procedure is called once per run.
18     ! It reads the pressure levels from "coefoz_LMDZ.nc".
19     ! We assume that, in "coefoz_LMDZ.nc", the pressure levels are in hPa
20     ! and strictly increasing.
21    
22     use netcdf95, only: nf95_open, nf95_close, nf95_get_coord
23     use netcdf, only: nf90_nowrite
24    
25     ! Variables local to the procedure:
26     integer ncid ! for NetCDF
27     integer n_plev ! number of pressure levels in the input data
28     integer k
29    
30     !---------------------------------------
31    
32     print *, "Call sequence information: press_coefoz"
33    
34     call nf95_open("coefoz_LMDZ.nc", nf90_nowrite, ncid)
35    
36     call nf95_get_coord(ncid, "plev", plev)
37     ! Convert from hPa to Pa because "regr_pr_av" and "regr_pr_int"
38     ! require so:
39     plev = plev * 100.
40     n_plev = size(plev)
41    
42     call nf95_close(ncid)
43    
44     ! Compute edges of pressure intervals:
45     allocate(press_in_edg(n_plev + 1))
46     press_in_edg(1) = 0.
47     ! We choose edges halfway in logarithm:
48     forall (k = 2:n_plev) press_in_edg(k) = sqrt(plev(k - 1) * plev(k))
49     press_in_edg(n_plev + 1) = huge(0.)
50     ! (infinity, but any value guaranteed to be greater than the
51     ! surface pressure would do)
52    
53     end subroutine press_coefoz
54    
55     end module press_coefoz_m

  ViewVC Help
Powered by ViewVC 1.1.21