/[lmdze]/trunk/phylmd/ozonecm.f
ViewVC logotype

Diff of /trunk/phylmd/ozonecm.f

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

trunk/libf/phylmd/ozonecm.f revision 21 by guez, Wed Feb 27 13:16:39 2008 UTC trunk/libf/phylmd/ozonecm.f90 revision 22 by guez, Fri Jul 31 15:18:47 2009 UTC
# Line 1  Line 1 
1  !  module ozonecm_m
2  ! $Header: /home/cvsroot/LMDZ4/libf/phylmd/ozonecm.F,v 1.3 2005/06/06 13:16:33 fairhead Exp $  
3  !    IMPLICIT NONE
4        SUBROUTINE ozonecm(rjour, rlat, paprs, o3)  
5        use dimens_m  contains
6        use dimphy  
7        use clesphys    function ozonecm(rjour, paprs)
8        use YOMCST  
9        IMPLICIT none      ! From phylmd/ozonecm.F, version 1.3 2005/06/06 13:16:33
10  C  
11  C The ozone climatology is based on an analytic formula which fits the      ! The ozone climatology is based on an analytic formula which fits the
12  C Krueger and Mintzner (1976) profile, as well as the variations with      ! Krueger and Mintzner (1976) profile, as well as the variations with
13  C altitude and latitude of the maximum ozone concentrations and the total      ! altitude and latitude of the maximum ozone concentrations and the total
14  C column ozone concentration of Keating and Young (1986). The analytic      ! column ozone concentration of Keating and Young (1986). The analytic
15  C formula have been established by J-F Royer (CRNM, Meteo France), who      ! formula have been established by J.-F. Royer (CRNM, Meteo France), who
16  C also provided us the code.      ! also provided us the code.
17  C  
18  C A. J. Krueger and R. A. Minzner, A Mid-Latitude Ozone Model for the      ! A. J. Krueger and R. A. Minzner, A Mid-Latitude Ozone Model for the
19  C 1976 U.S. Standard Atmosphere, J. Geophys. Res., 81, 4477, (1976).      ! 1976 U.S. Standard Atmosphere, J. Geophys. Res., 81, 4477, (1976).
20  C  
21  C Keating, G. M. and D. F. Young, 1985: Interim reference models for the      ! Keating, G. M. and D. F. Young, 1985: Interim reference models for the
22  C middle atmosphere, Handbook for MAP, vol. 16, 205-229.      ! middle atmosphere, Handbook for MAP, vol. 16, 205-229.
23  C  
24        use dimens_m, only: llm
25        real, intent(in):: rjour      USE dimphy, ONLY : klon
26        REAL, intent(in):: rlat(klon)      use numer_rec, only: assert, pi
27        real, intent(in):: paprs(klon,klev+1)      use phyetat0_m, only: rlat
28        REAL, intent(out):: o3(klon,klev)   ! ozone concentration in kg/kg  
29        REAL, INTENT (IN) :: rjour
30        REAL tozon      REAL, INTENT (IN) :: paprs(:, :) ! (klon, llm+1)
31        real pi, pl  
32        INTEGER i, k      REAL ozonecm(klon, llm)
33  C----------------------------------------------------------      ! "ozonecm(j, k)" is the column-density of ozone in cell "(j, k)", that is
34        REAL field(klon,klev+1)      ! between interface "k" and interface "k + 1", in kDU.
35        REAL ps  
36        PARAMETER (ps=101325.0)      ! Variables local to the procedure:
37        REAL an, unit, zo3q3  
38        SAVE an, unit, zo3q3      REAL tozon, pl
39        REAL mu,gms, zslat, zsint, zcost, z, ppm, qpm, a      INTEGER i, k
40        REAL asec, bsec, aprim, zo3a3  
41  C----------------------------------------------------------      REAL field(klon, llm+1)
42  c         data an /365.25/   (meteo)      ! "field(:, k)" is the column-density of ozone between interface
43        DATA an /360.00/      ! "k" and the top of the atmosphere (interface "llm + 1"), in kDU.
44        DATA unit /2.1415e-05/  
45        DATA zo3q3 /4.0e-08/      real, PARAMETER:: ps=101325.
46        REAL, parameter:: an = 360., unit = 2.1415E-5, zo3q3 = 4E-8
47        pi = 4.0 * ATAN(1.0)      REAL gms, zslat, zsint, zcost, z, ppm, qpm, a
48        DO k = 1, klev      REAL asec, bsec, aprim, zo3a3
49        DO i = 1, klon  
50        zslat = SIN(pi/180.*rlat(i))      !----------------------------------------------------------
51        zsint = SIN(2.*pi*(rjour+15.)/an)  
52        zcost = COS(2.*pi*(rjour+15.)/an)      call assert(shape(paprs) == (/klon, llm + 1/), "ozonecm")
53        z = 0.0531+zsint*(-0.001595+0.009443*zslat) +  
54       .  zcost*(-0.001344-0.00346*zslat) +      DO k = 1, llm
55       .  zslat**2*(.056222+zslat**2*(-.037609         DO i = 1, klon
56       . +.012248*zsint+.00521*zcost+.008890*zslat))            zslat = sin(pi / 180. * rlat(i))
57        zo3a3 = zo3q3/ps/2.            zsint = sin(2.*pi*(rjour+15.)/an)
58        z = z-zo3q3*ps            zcost = cos(2.*pi*(rjour+15.)/an)
59        gms = z            z = 0.0531 + zsint * (-0.001595+0.009443*zslat) &
60        mu = ABS(sin(pi/180.*rlat(i)))                 + zcost * (-0.001344-0.00346*zslat) &
61        ppm = 800.-(500.*zslat+150.*zcost)*zslat                 + zslat**2 * (.056222 + zslat**2 &
62        qpm = 1.74e-5-(7.5e-6*zslat+1.7e-6*zcost)*zslat                 * (-.037609+.012248*zsint+.00521*zcost+.008890*zslat))
63        bsec = 2650.+5000.*zslat**2            zo3a3 = zo3q3/ps/2.
64        a = 4.0*(bsec)**(3./2.)*(ppm)**(3./2.)*(1.0+(bsec/ps)**(3./2.))            z = z - zo3q3*ps
65        a = a/(bsec**(3./2.)+ppm**(3./2.))**2            gms = z
66        aprim = (2.666666*qpm*ppm-a*gms)/(1.0-a)            ppm = 800. - (500.*zslat+150.*zcost)*zslat
67        aprim = amax1(0.,aprim)            qpm = 1.74E-5 - (7.5E-6*zslat+1.7E-6*zcost)*zslat
68        asec = (gms-aprim)*(1.0+(bsec/ps)**(3./2.))            bsec = 2650. + 5000.*zslat**2
69        asec = AMAX1(0.0,asec)            a = 4.0*(bsec)**(3./2.)*(ppm)**(3./2.)*(1.0+(bsec/ps)**(3./2.))
70        aprim = gms-asec/(1.+(bsec/ps)**(3./2.))            a = a/(bsec**(3./2.)+ppm**(3./2.))**2
71        pl = paprs(i,k)            aprim = (2.666666*qpm*ppm-a*gms)/(1.0-a)
72        tozon = aprim/(1.+3.*(ppm/pl)**2)+asec/(1.+(bsec/pl)**(3./2.))            aprim = amax1(0., aprim)
73       .  + zo3a3*pl*pl            asec = (gms-aprim)*(1.0+(bsec/ps)**(3./2.))
74        tozon = tozon / 9.81  ! en kg/m**2            asec = amax1(0.0, asec)
75        tozon = tozon / unit  ! en kg/m**2  > u dobson (10e-5 m)            aprim = gms - asec/(1.+(bsec/ps)**(3./2.))
76        tozon = tozon / 1000. ! en cm            pl = paprs(i, k)
77        field(i,k) = tozon            tozon = aprim / (1. + 3. * (ppm / pl)**2) &
78        ENDDO                 + asec / (1. + (bsec / pl)**(3./2.)) + zo3a3 * pl * pl
79        ENDDO            ! Convert from Pa to kDU:
80  c            field(i, k) = tozon / 9.81 / unit / 1e3
81        DO i = 1, klon         END DO
82           field(i,klev+1) = 0.0      END DO
83        ENDDO  
84        DO k = 1, klev      field(:, llm+1) = 0.
85          DO i = 1, klon      forall (k = 1: llm) ozonecm(:, k) = field(:, k) - field(:, k+1)
86            o3(i,k) = field(i,k) - field(i,k+1)  
87            IF (.not. bug_ozone) then    END function ozonecm
88  c         convert o3 into kg/kg          
89              o3(i,k)=MAX(o3(i,k),1.0e-12)*RG/46.6968  end module ozonecm_m
      .               /(paprs(i,k)-paprs(i,k+1))  
           ENDIF  
         ENDDO  
       ENDDO  
 c  
       RETURN  
       END  

Legend:
Removed from v.21  
changed lines
  Added in v.22

  ViewVC Help
Powered by ViewVC 1.1.21