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

Diff of /trunk/phylmd/Mobidic/regr_pr_comb_coefoz.f

Parent Directory Parent Directory | Revision Log Revision Log | View Patch Patch

trunk/libf/phylmd/read_coefoz_m.f90 revision 3 by guez, Wed Feb 27 13:16:39 2008 UTC trunk/libf/phylmd/Mobidic/regr_pr_comb_coefoz.f90 revision 11 by guez, Thu Jun 5 12:43:08 2008 UTC
# Line 1  Line 1 
1  module read_coefoz_m  module regr_pr_comb_coefoz_m
2    
3    ! This module is clean: no C preprocessor directive, no include line.    ! This module is clean: no C preprocessor directive, no include line.
4    
# Line 7  module read_coefoz_m Line 7  module read_coefoz_m
7    
8    implicit none    implicit none
9    
10    real, save:: c_Mob(klon, llm, 12)    ! The five module variables declared here are on the "physics" grid.
11      ! The value of each variable for index "(i, k)" is at longitude
12      ! "rlon(i)", latitude "rlat(i)" and middle of layer "k".
13    
14      real, save:: c_Mob(klon, llm)
15    ! (sum of Mobidic terms in the net mass production rate of ozone    ! (sum of Mobidic terms in the net mass production rate of ozone
16    ! by chemistry, per unit mass of air, in s-1)    ! by chemistry, per unit mass of air, in s-1)
   ! (On the "physics" grid.  
   ! Third dimension is the number of the month in the year.  
   ! "c_Mob(i, k, month)" is at longitude "rlon(i)", latitude "rlat(i)",  
   ! middle of layer "k".)  
17    
18    real, save:: a2(klon, llm, 12)    real, save:: a2(klon, llm)
19    ! (derivative of mass production rate of ozone per unit mass of    ! (derivative of mass production rate of ozone per unit mass of
20    ! air with respect to ozone mass fraction, in s-1)    ! air with respect to ozone mass fraction, in s-1)
   ! (On the "physics" grid.  
   ! Third dimension is the number of the month in the year.  
   ! "a2(i, k, month)" is at longitude "rlon(i)", latitude "rlat(i)",  
   ! middle of layer "k".)  
21    
22    real, save:: a4_mass(klon, llm, 12)    real, save:: a4_mass(klon, llm)
23    ! (derivative of mass production rate of ozone per unit mass of    ! (derivative of mass production rate of ozone per unit mass of
24    ! air with respect to temperature, in s-1 K-1)    ! air with respect to temperature, in s-1 K-1)
   ! (On the "physics" grid.  
   ! Third dimension is the number of the month in the year.  
   ! "a4_mass(i, k, month)" is at longitude "rlon(i)", latitude "rlat(i)",  
   ! middle of layer "k".)  
25    
26    real, save:: a6_mass(klon, llm, 12)    real, save:: a6_mass(klon, llm)
27    ! (derivative of mass production rate of ozone per unit mass of    ! (derivative of mass production rate of ozone per unit mass of
28    ! air with respect to mass column-density of ozone above, in m2 s-1 kg-1)    ! air with respect to mass column-density of ozone above, in m2 s-1 kg-1)
   ! (On the "physics" grid.  
   ! Third dimension is the number of the month in the year.  
   ! "a6_mass(i, k, month)" is at longitude "rlon(i)", latitude "rlat(i)",  
   ! middle of layer "k".)  
29    
30    real, save:: r_het_interm(klon, llm, 12)    real, save:: r_het_interm(klon, llm)
31    ! (net mass production rate by heterogeneous chemistry, per unit    ! (net mass production rate by heterogeneous chemistry, per unit
32    ! mass of ozone, corrected for chlorine content and latitude, but    ! mass of ozone, corrected for chlorine content and latitude, but
33    ! not for temperature and sun direction, in s-1)    ! not for temperature and sun direction, in s-1)
   ! (On the "physics" grid.  
   ! Third dimension is the number of the month in the year.  
   ! "r_het_interm(i, k, month)" is at longitude "rlon(i)", latitude "rlat(i)",  
   ! middle of layer "k".)  
34    
35    private klon, llm    private klon, llm
36    
37  contains  contains
38    
39    subroutine read_coefoz    subroutine regr_pr_comb_coefoz(julien)
40    
41      ! This subroutine reads from a file all eight parameters for ozone      ! "regr_pr_comb_coefoz" stands for "regrid pressure combine
42      ! chemistry, at all months.      ! coefficients ozone".
     ! The parameters are packed to the "physics" grid.  
     ! Finally, the eight parameters are combined to define the five  
     ! module variables.  
43    
44      use netcdf95, only: nf95_open, nf95_close, nf90_nowrite      ! This subroutine :
45      use o3_Mob_ph_m, only: o3_Mob_ph      ! -- reads from a file all eight coefficients for ozone chemistry,
46        !    at the current day ;
47        ! -- regrids the coefficients in pressure to the LMDZ vertical grid ;
48        ! -- packs the coefficients to the "physics" horizontal grid ;
49        ! -- combines the eight coefficients to define the five module variables.
50    
51        ! We assume that, in "coefoz_LMDZ.nc", the pressure levels are in hPa
52        ! and strictly increasing.
53    
54        use netcdf95, only: nf95_open, nf95_close, nf95_get_coord
55        use netcdf, only: nf90_nowrite
56        use regr_pr_coefoz, only: regr_pr_av_coefoz, regr_pr_int_coefoz
57      use phyetat0_m, only: rlat      use phyetat0_m, only: rlat
58    
59        integer, intent(in):: julien ! jour julien, 1 <= julien <= 360
60    
61      ! Variables local to the procedure:      ! Variables local to the procedure:
62      integer ncid ! for NetCDF      integer ncid ! for NetCDF
63    
64      real a6(klon, llm, 12)      real, pointer:: plev(:)
65      ! (derivative of P_net_Mob with respect to column-density of ozone      ! (pressure level of input data, converted to Pa, in strictly
66        ! increasing order)
67    
68        integer n_plev ! number of pressure levels in the input data
69    
70        real, allocatable:: press_in_edg(:)
71        ! (edges of pressure intervals for input data, in Pa, in strictly
72        ! increasing order)
73    
74        real coefoz(klon, llm)
75        ! (temporary storage for an ozone coefficient)
76        ! (On the "physics" grid.
77        ! "coefoz(i, k)" is at longitude "rlon(i)", latitude "rlat(i)",
78        ! middle of layer "k".)
79    
80        real a6(klon, llm)
81        ! (derivative of "P_net_Mob" with respect to column-density of ozone
82      ! above, in cm2 s-1)      ! above, in cm2 s-1)
83      ! (On the "physics" grid.      ! (On the "physics" grid.
84      ! Third dimension is the number of the month in the year.      ! "a6(i, k)" is at longitude "rlon(i)", latitude "rlat(i)",
     ! "a6(i, k, month)" is at longitude "rlon(i)", latitude "rlat(i)",  
85      ! middle of layer "k".)      ! middle of layer "k".)
86    
87      real, parameter:: amu = 1.6605402e-27 ! atomic mass unit, in kg      real, parameter:: amu = 1.6605402e-27 ! atomic mass unit, in kg
# Line 80  contains Line 89  contains
89      real, parameter:: Clx = 3.8e-9      real, parameter:: Clx = 3.8e-9
90      ! (total chlorine content in the upper stratosphere)      ! (total chlorine content in the upper stratosphere)
91    
92      integer k, month      integer k
93    
94      !------------------------------------      !------------------------------------
95    
96      print *, "Call sequence information: read_coefoz"      print *, "Call sequence information: regr_pr_comb_coefoz"
97    
98      call nf95_open("coefoz_LMDZ.nc", nf90_nowrite, ncid)      call nf95_open("coefoz_LMDZ.nc", nf90_nowrite, ncid)
99    
100      a2 = o3_Mob_ph(ncid, "a2")      call nf95_get_coord(ncid, "plev", plev)
101      a4_mass = o3_Mob_ph(ncid, "a4") * 48. / 29.      ! Convert from hPa to Pa because "regr_pr_av" and "regr_pr_int" require so:
102      a6 = o3_Mob_ph(ncid, "a6")      plev = plev * 100.
103        n_plev = size(plev)
104    
105        ! Compute edges of pressure intervals:
106        allocate(press_in_edg(n_plev + 1))
107        press_in_edg(1) = 0.
108        ! We choose edges halfway in logarithm:
109        forall (k = 2:n_plev) press_in_edg(k) = sqrt(plev(k - 1) * plev(k))
110        press_in_edg(n_plev + 1) = huge(0.)
111        ! (infinity, but any value guaranteed to be greater than the
112        ! surface pressure would do)
113    
114        call regr_pr_av_coefoz(ncid, "a2", julien, press_in_edg, a2)
115    
116        call regr_pr_av_coefoz(ncid, "a4", julien, press_in_edg, a4_mass)
117        a4_mass = a4_mass * 48. / 29.
118    
119        call regr_pr_av_coefoz(ncid, "a6", julien, press_in_edg, a6)
120    
121      ! Compute "a6_mass" avoiding underflow, do not divide by 1e4      ! Compute "a6_mass" avoiding underflow, do not divide by 1e4
122      ! before dividing by molecular mass:      ! before dividing by molecular mass:
123      a6_mass = a6 / (1e4 * 29. * amu)      a6_mass = a6 / (1e4 * 29. * amu)
124      ! (factor 1e4: conversion from cm2 to m2)      ! (factor 1e4: conversion from cm2 to m2)
125    
126      c_Mob = 48. / 29. * (o3_Mob_ph(ncid, "P_net_Mob") &      ! Combine coefficients to get "c_Mob":
127           - a2 * o3_Mob_ph(ncid, "r_Mob") - a6 * o3_Mob_ph(ncid, "Sigma_Mob")) &      ! (We use as few local variables as possible, in order to spare
128           - a4_mass * o3_Mob_ph(ncid, "temp_Mob")      ! main memory.)
129    
130        call regr_pr_av_coefoz(ncid, "P_net_Mob", julien, press_in_edg, c_Mob)
131    
132        call regr_pr_av_coefoz(ncid, "r_Mob", julien, press_in_edg, coefoz)
133        c_mob = c_mob - a2 * coefoz
134    
135        call regr_pr_int_coefoz(ncid, "Sigma_Mob", julien, plev, top_value=0., &
136             v3=coefoz)
137        c_mob = (c_mob - a6 * coefoz) * 48. / 29.
138    
139        call regr_pr_av_coefoz(ncid, "temp_Mob", julien, press_in_edg, coefoz)
140        c_mob = c_mob - a4_mass * coefoz
141    
142      r_het_interm = o3_Mob_ph(ncid, "R_Het") * (Clx / 3.8e-9)**2      call regr_pr_av_coefoz(ncid, "R_Het", julien, press_in_edg, r_het_interm)
143      ! Heterogeneous chemistry is only at high latitudes:      ! Heterogeneous chemistry is only at high latitudes:
144      forall (k = 1: llm, month = 1: 12)      forall (k = 1: llm)
145         where (abs(rlat) <= 45.) r_het_interm(:, k, month) = 0.         where (abs(rlat) <= 45.) r_het_interm(:, k) = 0.
146      end forall      end forall
147        r_het_interm = r_het_interm * (Clx / 3.8e-9)**2
148    
149        deallocate(plev) ! pointer
150      call nf95_close(ncid)      call nf95_close(ncid)
151    
152    end subroutine read_coefoz    end subroutine regr_pr_comb_coefoz
153    
154  end module read_coefoz_m  end module regr_pr_comb_coefoz_m

Legend:
Removed from v.3  
changed lines
  Added in v.11

  ViewVC Help
Powered by ViewVC 1.1.21