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

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

Parent Directory Parent Directory | Revision Log Revision Log


Revision 20 - (hide annotations)
Wed Oct 15 16:19:57 2008 UTC (15 years, 7 months ago) by guez
File size: 4204 byte(s)
Deleted argument "presnivs" of "physiq", "ini_histhf", "ini_histhf3d",
"ini_histday", "ini_histins", "ini_histrac", "phytrac". Access it from
"comvert" instead.

Replaced calls to NetCDF Fortran 77 interface by calls to Fortran 90
interface or to NetCDF95.

Procedure "gr_phy_write_3d" now works with a variable of arbitrary
size in the second dimension.

Annotated use statements with "only" clause.

Replaced calls to NetCDF interface version 2 by calls to Fortran 90
interface in "guide.f90" and "read_reanalyse.f".

In "write_histrac", replaced calls to "gr_fi_ecrit" by calls to
"gr_phy_write_2d" and "gr_phy_write_3d".

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 guez 19 ! Ozone mole fraction from "coefoz_LMDZ.nc" is a 2D latitude --
13     ! pressure variable.
14     ! The target horizontal LMDZ grid is the "scalar" grid: "rlonv", "rlatu".
15     ! The target vertical LMDZ grid is the grid of layer boundaries.
16     ! We assume that the input variable is already on the LMDZ "rlatu"
17     ! latitude grid.
18     ! The input variable does not depend on longitude, but the
19     ! pressure at LMDZ layers does.
20     ! Therefore, the values on the LMDZ grid do depend on longitude.
21 guez 8 ! Regridding is by averaging, assuming a step function.
22 guez 19 ! We assume that, in the input file, the pressure levels are in
23     ! hPa and strictly increasing.
24 guez 8
25     use conf_gcm_m, only: dayref
26     use dimens_m, only: iim, jjm, llm
27     use netcdf95, only: nf95_open, nf95_close, nf95_inq_varid, handle_err, &
28     nf95_get_coord
29     use netcdf, only: nf90_nowrite, nf90_get_var
30 guez 13 use numer_rec, only: assert
31 guez 19 use grid_change, only: dyn_phy
32     use regr1_step_av_m, only: regr1_step_av
33     use pressure_var, only: p3d
34 guez 8
35     real, intent(out):: o3_mob_regr(:, :, :) ! (iim + 1, jjm + 1, llm)
36     ! (ozone mole fraction from Mobidic adapted to the LMDZ grid)
37     ! ("o3_mob_regr(i, j, l)" is at longitude "rlonv(i)", latitude
38     ! "rlatu(j)" and pressure level "pls(i, j, l)")
39    
40     ! Variables local to the procedure:
41    
42     real, pointer:: plev(:)
43     ! (pressure levels of Mobidic data, in Pa, in strictly increasing order)
44    
45     real, allocatable:: press_in_edg(:)
46     ! (edges of pressure intervals for Mobidic data, in Pa, in strictly
47     ! increasing order)
48    
49     integer ncid, varid, ncerr ! for NetCDF
50     integer n_plev ! number of pressure levels in Mobidic data
51 guez 19 integer i, j
52 guez 8
53     real, allocatable:: r_mob(:, :)! (jjm + 1, n_plev)
54     ! (ozone mole fraction from Mobidic at day "dayref")
55     ! (r_mob(j, k) is at latitude "rlatu(j)" and pressure level "plev(k)".)
56    
57     !------------------------------------------------------------
58    
59 guez 10 print *, "Call sequence information: regr_pr_o3"
60 guez 8 call assert(shape(o3_mob_regr) == (/iim + 1, jjm + 1, llm/), "regr_pr_o3")
61    
62     call nf95_open("coefoz_LMDZ.nc", nf90_nowrite, ncid)
63    
64     call nf95_get_coord(ncid, "plev", plev)
65     ! Convert from hPa to Pa because "regr_pr_av" requires so:
66     plev = plev * 100.
67     n_plev = size(plev)
68    
69     ! Compute edges of pressure intervals:
70     allocate(press_in_edg(n_plev + 1))
71     press_in_edg(1) = 0.
72     ! We choose edges halfway in logarithm:
73     forall (j = 2:n_plev) press_in_edg(j) = sqrt(plev(j - 1) * plev(j))
74     press_in_edg(n_plev + 1) = huge(0.)
75     ! (infinity, but any value guaranteed to be greater than the
76     ! surface pressure would do)
77    
78     deallocate(plev) ! pointer
79    
80     call nf95_inq_varid(ncid, "r_Mob", varid)
81     allocate(r_mob(jjm + 1, n_plev))
82    
83     ! Get data at the right day from the input file:
84     ncerr = nf90_get_var(ncid, varid, r_mob, start=(/1, 1, dayref/))
85     call handle_err("nf90_get_var r_Mob", ncerr)
86     ! Latitudes are in increasing order in the input file while
87     ! "rlatu" is in decreasing order so we need to invert order:
88     r_mob = r_mob(jjm+1:1:-1, :)
89    
90     call nf95_close(ncid)
91    
92 guez 20 ! Regrid in pressure by averaging a step function of pressure:
93 guez 19 do j = 1, jjm + 1
94     do i = 1, iim
95     if (dyn_phy(i, j)) then
96     o3_mob_regr(i, j, llm:1:-1) &
97     = regr1_step_av(r_mob(j, :), press_in_edg, &
98     p3d(i, j, llm+1:1:-1))
99     ! (invert order of indices because "p3d" is decreasing)
100     end if
101     end do
102     end do
103 guez 8
104 guez 19 ! Duplicate pole values on all longitudes:
105     o3_mob_regr(2:, 1, :) = spread(o3_mob_regr(1, 1, :), dim=1, ncopies=iim)
106     o3_mob_regr(2:, jjm + 1, :) &
107     = spread(o3_mob_regr(1, jjm + 1, :), dim=1, ncopies=iim)
108    
109     ! Duplicate first longitude to last longitude:
110     o3_mob_regr(iim + 1, 2:jjm, :) = o3_mob_regr(1, 2:jjm, :)
111    
112 guez 8 end subroutine regr_pr_o3
113    
114     end module regr_pr_o3_m

  ViewVC Help
Powered by ViewVC 1.1.21