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

Diff of /trunk/Sources/phylmd/ozonecm.f

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

trunk/libf/phylmd/ozonecm.f revision 3 by guez, Wed Feb 27 13:16:39 2008 UTC trunk/Sources/phylmd/ozonecm.f revision 134 by guez, Wed Apr 29 15:47:56 2015 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 nr_util, 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  
31        real pi, pl      REAL, INTENT (IN) :: paprs(:, :) ! (klon, llm+1)
32        INTEGER i, k      ! pression pour chaque inter-couche, en Pa
33  C----------------------------------------------------------  
34        REAL field(klon,klev+1)      REAL ozonecm(klon, llm)
35        REAL ps      ! "ozonecm(j, k)" is the column-density of ozone in cell "(j, k)", that is
36        PARAMETER (ps=101325.0)      ! between interface "k" and interface "k + 1", in kDU.
37        REAL an, unit, zo3q3  
38        SAVE an, unit, zo3q3      ! Local:
39        REAL mu,gms, zslat, zsint, zcost, z, ppm, qpm, a  
40        REAL asec, bsec, aprim, zo3a3      REAL tozon ! equivalent pressure of ozone above interface "k", in Pa
41  C----------------------------------------------------------      INTEGER i, k
42  c         data an /365.25/   (meteo)  
43        DATA an /360.00/      REAL field(llm+1)
44        DATA unit /2.1415e-05/      ! "field(k)" is the column-density of ozone between interface
45        DATA zo3q3 /4.0e-08/      ! "k" and the top of the atmosphere (interface "llm + 1"), in kDU.
46    
47        pi = 4.0 * ATAN(1.0)      real, PARAMETER:: ps = 101325.
48        DO k = 1, klev      REAL, parameter:: an = 360., zo3q3 = 4E-8
49        DO i = 1, klon      real, parameter:: zo3a3 = zo3q3 / ps / 2.
50        zslat = SIN(pi/180.*rlat(i))      REAL, parameter:: dobson_unit = 2.1415E-5 ! in kg m-2
51        zsint = SIN(2.*pi*(rjour+15.)/an)      REAL gms, slat, slat2, sint, cost, ppm, a
52        zcost = COS(2.*pi*(rjour+15.)/an)      REAL asec, bsec, aprim
53        z = 0.0531+zsint*(-0.001595+0.009443*zslat) +  
54       .  zcost*(-0.001344-0.00346*zslat) +      !----------------------------------------------------------
55       .  zslat**2*(.056222+zslat**2*(-.037609  
56       . +.012248*zsint+.00521*zcost+.008890*zslat))      call assert(shape(paprs) == (/klon, llm + 1/), "ozonecm")
57        zo3a3 = zo3q3/ps/2.  
58        z = z-zo3q3*ps      sint = sin(2 * pi * (rjour + 15.) / an)
59        gms = z      cost = cos(2 * pi * (rjour + 15.) / an)
60        mu = ABS(sin(pi/180.*rlat(i)))      field(llm + 1) = 0.
61        ppm = 800.-(500.*zslat+150.*zcost)*zslat  
62        qpm = 1.74e-5-(7.5e-6*zslat+1.7e-6*zcost)*zslat      DO i = 1, klon
63        bsec = 2650.+5000.*zslat**2         slat = sin(pi / 180. * rlat(i))
64        a = 4.0*(bsec)**(3./2.)*(ppm)**(3./2.)*(1.0+(bsec/ps)**(3./2.))         slat2 = slat * slat
65        a = a/(bsec**(3./2.)+ppm**(3./2.))**2         gms = 0.0531 + sint * (- 0.001595 + 0.009443 * slat) + cost &
66        aprim = (2.666666*qpm*ppm-a*gms)/(1.0-a)              * (- 0.001344 - 0.00346 * slat) + slat2 * (0.056222 + slat2 &
67        aprim = amax1(0.,aprim)              * (- 0.037609 + 0.012248 * sint + 0.00521 * cost + 0.00889 &
68        asec = (gms-aprim)*(1.0+(bsec/ps)**(3./2.))              * slat)) - zo3q3 * ps
69        asec = AMAX1(0.0,asec)         ppm = 800. - 500. * slat2 - 150. * cost * slat
70        aprim = gms-asec/(1.+(bsec/ps)**(3./2.))         bsec = 2650. + 5000. * slat2
71        pl = paprs(i,k)         a = 4. * bsec**1.5 * ppm**1.5 * (1. + (bsec / ps)**1.5) &
72        tozon = aprim/(1.+3.*(ppm/pl)**2)+asec/(1.+(bsec/pl)**(3./2.))              / (bsec**1.5 + ppm**1.5)**2
73       .  + zo3a3*pl*pl         aprim = max(0., (2.666666 * (1.74E-5 - 7.5E-6 * slat2 &
74        tozon = tozon / 9.81  ! en kg/m**2              - 1.7E-6 * cost * slat) * ppm - a * gms) / (1. - a))
75        tozon = tozon / unit  ! en kg/m**2  > u dobson (10e-5 m)         asec = max(0., (gms - aprim) * (1. + (bsec / ps)**1.5))
76        tozon = tozon / 1000. ! en cm         aprim = gms - asec / (1. + (bsec / ps)**1.5)
77        field(i,k) = tozon  
78        ENDDO         DO k = 1, llm
79        ENDDO            tozon = aprim / (1. + 3. * (ppm / paprs(i, k))**2) &
80  c                 + asec / (1. + (bsec / paprs(i, k))**1.5) &
81        DO i = 1, klon                 + zo3a3 * paprs(i, k) * paprs(i, k)
82           field(i,klev+1) = 0.0            ! Convert from Pa to kDU:
83        ENDDO            field(k) = tozon / 9.81 / dobson_unit / 1e3
84        DO k = 1, klev         END DO
85          DO i = 1, klon  
86            o3(i,k) = field(i,k) - field(i,k+1)         forall (k = 1: llm) ozonecm(i, k) = field(k) - field(k + 1)
87            IF (.not. bug_ozone) then      END DO
88  c         convert o3 into kg/kg          
89              o3(i,k)=MAX(o3(i,k),1.0e-12)*RG/46.6968      ozonecm = max(ozonecm, 1e-12)
90       .               /(paprs(i,k)-paprs(i,k+1))  
91            ENDIF    END function ozonecm
92          ENDDO  
93        ENDDO  end module ozonecm_m
 c  
       RETURN  
       END  

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

  ViewVC Help
Powered by ViewVC 1.1.21