source: trunk/LMDZ.GENERIC/libf/phystd/aeropacity.F90 @ 135

Last change on this file since 135 was 135, checked in by aslmd, 14 years ago

CHANGEMENT ARBORESCENCE ETAPE 2 -- NON COMPLET

File size: 5.7 KB
Line 
1      subroutine aeropacity(ngrid,nlayer,nq,pplev,pq, &
2         aerosol,reffrad,QREFvis3d,QREFir3d,tau_col)
3
4       use radinc_h, only : naerkind
5                 
6       implicit none
7
8!==================================================================
9!     
10!     Purpose
11!     -------
12!     Compute aerosol optical depth in each gridbox.
13!     
14!     Authors
15!     -------
16!     F. Forget
17!     F. Montmessin (water ice scheme)
18!     update J.-B. Madeleine (2008)
19!     dust removal, simplification by Robin Wordsworth (2009)
20!
21!     Input
22!     -----
23!     ngrid             Number of horizontal gridpoints
24!     nlayer            Number of layers
25!     nq                Number of tracers
26!     pplev             Pressure (Pa) at each layer boundary
27!     pq                Aerosol mixing ratio
28!     reffrad(ngrid,nlayer,naerkind)         Aerosol effective radius
29!     QREFvis3d(ngridmx,nlayermx,naerkind) \ 3d extinction coefficients
30!     QREFir3d(ngridmx,nlayermx,naerkind)  / at reference wavelengths
31!
32!     Output
33!     ------
34!     aerosol            Aerosol optical depth in layer l, grid point ig
35!     tau_col            Total column optical depth at grid point ig
36!
37!=======================================================================
38
39#include "dimensions.h"
40#include "dimphys.h"
41#include "callkeys.h"
42#include "comcstfi.h"
43#include "comgeomfi.h"
44#include "tracer.h"
45
46      INTEGER ngrid,nlayer,nq
47
48      REAL pplev(ngrid,nlayer+1)
49      REAL pq(ngrid,nlayer,nq)
50      REAL aerosol(ngrid,nlayer,naerkind)
51      REAL reffrad(ngrid,nlayer,naerkind)
52      REAL QREFvis3d(ngridmx,nlayermx,naerkind)
53      REAL QREFir3d(ngridmx,nlayermx,naerkind)
54
55      REAL tauref(ngrid), tau_col(ngrid)
56
57      INTEGER l,ig,iq,iaer
58
59      LOGICAL firstcall
60      DATA firstcall/.true./
61      SAVE firstcall
62
63      REAL CBRT
64      EXTERNAL CBRT
65
66      INTEGER,SAVE :: i_co2ice=0      ! co2 ice
67      INTEGER,SAVE :: i_h2oice=0      ! water ice
68      CHARACTER(LEN=20) :: tracername ! to temporarily store text
69
70
71! identify tracers
72      IF (firstcall) THEN
73
74         ! are these tests of any real use?
75        IF(ngrid.NE.ngridmx) THEN
76            PRINT*,'STOP in aeropacity!'
77            PRINT*,'problem of dimensions:'
78            PRINT*,'ngrid  =',ngrid
79            PRINT*,'ngridmx  =',ngridmx
80            STOP
81        ENDIF
82
83        if (nq.gt.nqmx) then
84           write(*,*) 'STOP in aeropacity: (nq .gt. nqmx)!'
85           write(*,*) 'nq=',nq,' nqmx=',nqmx
86           stop
87        endif
88
89
90        do iq=1,nqmx
91          tracername=noms(iq)
92          if (tracername.eq."co2_ice") then
93            i_co2ice=iq
94          endif
95          if (tracername.eq."h2o_ice") then
96            i_h2oice=iq
97          endif
98        enddo
99       
100        write(*,*) "aeropacity: i_co2ice=",i_co2ice
101        write(*,*) "            i_h2oice=",i_h2oice
102       
103        !if((i_h2oice.lt.1) .and. (naerkind.gt.1))then
104        !   print*,'naerkind > 1 but no h2o ice, this will not work.'
105        !   stop
106        !endif
107
108        firstcall=.false.
109      ENDIF ! of IF (firstcall)
110
111
112      aerosol(:,:,:)=0.0
113
114      DO iaer = 1, naerkind ! Loop on aerosol kind (NOT tracer)
115!     ---------------------------------------------------------
116        aerkind: SELECT CASE (iaer)
117!==================================================================
118        CASE(1) aerkind                          ! CO2 ice (iaer=1)
119!==================================================================
120
121!       1. Initialization
122        CALL zerophys(ngrid*nlayer,aerosol(1,1,iaer))
123
124!       2. Opacity calculation
125        if (aerofixed.or.(i_co2ice.eq.0)) then !  CO2 ice cloud prescribed
126           do ig=1, ngrid
127              do l=1,nlayer
128                 aerosol(ig,l,iaer) =1.e-9
129              end do
130              aerosol(ig,14,iaer)=1.e-9 ! single cloud layer option
131           end do
132        else
133           DO ig=1, ngrid
134              DO l=1,nlayer
135                 aerosol(ig,l,iaer) =                   &
136                      (  0.75 * QREFvis3d(ig,l,iaer) /  &
137                      ( rho_co2 * reffrad(ig,l,iaer) )  ) * &
138                      ( pq(ig,l,i_co2ice) + 1.E-9 ) *   &
139                      ( pplev(ig,l) - pplev(ig,l+1) ) / g
140              ENDDO
141           ENDDO
142        end if
143
144!==================================================================
145        CASE(2) aerkind               ! Water ice crystals (iaer=2)
146!==================================================================
147
148
149!       1. Initialization
150        CALL zerophys(ngrid*nlayer,aerosol(1,1,iaer))
151
152!       2. Opacity calculation
153!     if (1.eq.1) then
154     if (aerofixed.or.(i_h2oice.eq.0)) then
155!        print*,'No H2O clouds in the radiative transfer!'
156        do ig=1, ngrid
157           do l=1,nlayer
158              aerosol(ig,l,iaer) =1.e-9
159           end do
160           aerosol(ig,5,iaer)=1.e-9 ! single cloud layer option
161        end do
162     else
163        DO ig=1, ngrid
164          DO l=1,nlayer
165             aerosol(ig,l,iaer) =                   &
166            (  0.75 * QREFvis3d(ig,l,iaer) /        &
167              ( rho_ice * reffrad(ig,l,iaer) )  ) * &
168              ( pq(ig,l,i_h2oice) + 1.E-9 ) *   &
169              ( pplev(ig,l) - pplev(ig,l+1) ) / g
170          ENDDO
171        ENDDO
172
173     end if
174
175
176!==================================================================
177        END SELECT aerkind
178
179      ENDDO ! iaer (loop on aerosol kind)
180
181
182! --------------------------------------------------------------------------
183! Column integrated visible optical depth in each point (used for diagnostic)
184
185      call zerophys(ngrid,tau_col)
186      do iaer = 1, naerkind
187         do l=1,nlayer
188            do ig=1,ngrid
189               tau_col(ig) = tau_col(ig) + aerosol(ig,l,iaer)
190            end do
191         end do
192      end do
193
194      return
195    end subroutine aeropacity
196     
Note: See TracBrowser for help on using the repository browser.