/[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 22 - (hide annotations)
Fri Jul 31 15:18:47 2009 UTC (14 years, 10 months ago) by guez
File size: 4242 byte(s)
Superficial modifications
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 guez 22 nf95_gw_var
29 guez 8 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 guez 22 call nf95_inq_varid(ncid, "plev", varid)
65     call nf95_gw_var(ncid, varid, plev)
66 guez 8 ! Convert from hPa to Pa because "regr_pr_av" requires so:
67     plev = plev * 100.
68     n_plev = size(plev)
69    
70     ! Compute edges of pressure intervals:
71     allocate(press_in_edg(n_plev + 1))
72     press_in_edg(1) = 0.
73     ! We choose edges halfway in logarithm:
74     forall (j = 2:n_plev) press_in_edg(j) = sqrt(plev(j - 1) * plev(j))
75     press_in_edg(n_plev + 1) = huge(0.)
76     ! (infinity, but any value guaranteed to be greater than the
77     ! surface pressure would do)
78    
79     deallocate(plev) ! pointer
80    
81     call nf95_inq_varid(ncid, "r_Mob", varid)
82     allocate(r_mob(jjm + 1, n_plev))
83    
84     ! Get data at the right day from the input file:
85     ncerr = nf90_get_var(ncid, varid, r_mob, start=(/1, 1, dayref/))
86     call handle_err("nf90_get_var r_Mob", ncerr)
87     ! Latitudes are in increasing order in the input file while
88     ! "rlatu" is in decreasing order so we need to invert order:
89     r_mob = r_mob(jjm+1:1:-1, :)
90    
91     call nf95_close(ncid)
92    
93 guez 20 ! Regrid in pressure by averaging a step function of pressure:
94 guez 19 do j = 1, jjm + 1
95     do i = 1, iim
96     if (dyn_phy(i, j)) then
97     o3_mob_regr(i, j, llm:1:-1) &
98     = regr1_step_av(r_mob(j, :), press_in_edg, &
99     p3d(i, j, llm+1:1:-1))
100     ! (invert order of indices because "p3d" is decreasing)
101     end if
102     end do
103     end do
104 guez 8
105 guez 19 ! Duplicate pole values on all longitudes:
106     o3_mob_regr(2:, 1, :) = spread(o3_mob_regr(1, 1, :), dim=1, ncopies=iim)
107     o3_mob_regr(2:, jjm + 1, :) &
108     = spread(o3_mob_regr(1, jjm + 1, :), dim=1, ncopies=iim)
109    
110     ! Duplicate first longitude to last longitude:
111     o3_mob_regr(iim + 1, 2:jjm, :) = o3_mob_regr(1, 2:jjm, :)
112    
113 guez 8 end subroutine regr_pr_o3
114    
115     end module regr_pr_o3_m

  ViewVC Help
Powered by ViewVC 1.1.21