source: trunk/LMDZ.PLUTO.old/libf/phypluto/optcv.F90 @ 3436

Last change on this file since 3436 was 3175, checked in by emillour, 11 months ago

Pluto PCM:
Add the old Pluto LMDZ for reference (required prior step to making
an LMDZ.PLUTO using the same framework as the other physics packages).
TB+EM

File size: 8.5 KB
Line 
1SUBROUTINE OPTCV(DTAUV,TAUV,TAUCUMV,PLEV,  &
2      QXVAER,QSVAER,GVAER,WBARV,COSBV,      &
3      TAURAY,TAUAERO,TMID,PMID,TAUGSURF,QVAR)
4
5      use radinc_h
6      use radcommon_h, only: gasv, tlimit, wrefVAR, Cmk, tgasref, &
7          pfgasref
8
9      implicit none
10
11!==================================================================
12!     
13!     Purpose
14!     -------
15!     Calculates shortwave optical constants at each level.
16!     
17!     Authors
18!     -------
19!     Adapted from the NASA Ames code by R. Wordsworth (2009)
20!     
21!==================================================================
22
23
24
25!     THIS SUBROUTINE SETS THE OPTICAL CONSTANTS IN THE VISUAL 
26!     IT CALCUALTES FOR EACH LAYER, FOR EACH SPECRAL INTERVAL IN THE VISUAL
27!     LAYER: WBAR, DTAU, COSBAR
28!     LEVEL: TAU
29!     
30!     TAUV(L,NW,NG) is the cumulative optical depth at the top of radiation code
31!     layer L. NW is spectral wavelength interval, ng the Gauss point index.
32!     
33!     TLEV(L) - Temperature at the layer boundary
34!     PLEV(L) - Pressure at the layer boundary (i.e. level)
35!     GASV(NT,NPS,NW,NG) - Visual CO2 k-coefficients
36!     
37!----------------------------------------------------------------------C
38
39#include "callkeys.h"
40#include "comcstfi.h"
41
42      real*8 DTAUV(L_NLAYRAD,L_NSPECTV,L_NGAUSS)
43      real*8 DTAUKV(L_LEVELS,L_NSPECTV,L_NGAUSS)
44      real*8 TAUV(L_NLEVRAD,L_NSPECTV,L_NGAUSS)
45      real*8 TAUCUMV(L_LEVELS,L_NSPECTV,L_NGAUSS)
46      real*8 PLEV(L_LEVELS)
47      real*8 TMID(L_LEVELS), PMID(L_LEVELS)
48      real*8 COSBV(L_NLAYRAD,L_NSPECTV,L_NGAUSS)
49      real*8 WBARV(L_NLAYRAD,L_NSPECTV,L_NGAUSS)
50
51!     For aerosols
52      real*8 QXVAER(L_LEVELS,L_NSPECTV,NAERKIND)
53      real*8 QSVAER(L_LEVELS,L_NSPECTV,NAERKIND)
54      real*8 GVAER(L_LEVELS,L_NSPECTV,NAERKIND)
55      real*8 TAUAERO(L_LEVELS,NAERKIND)
56      real*8 TAUAEROLK(L_LEVELS,L_NSPECTV,NAERKIND)
57      real*8 TAEROS(L_LEVELS,L_NSPECTV,NAERKIND)
58
59      integer L, NW, NG, K, NG1(L_NSPECTV), LK, IAER
60      integer MT(L_LEVELS), MP(L_LEVELS), NP(L_LEVELS)
61      real*8  ANS, TAUGAS
62      real*8 TAURAY(L_NSPECTV)
63      real*8  TRAY(L_LEVELS,L_NSPECTV)
64      real*8  DPR(L_LEVELS), U(L_LEVELS)
65      real*8  LCOEF(4), LKCOEF(L_LEVELS,4)
66
67      real*8 taugsurf(L_NSPECTV,L_NGAUSS-1), TRAYAER
68
69
70
71
72
73!     mixing ratio variables
74      real*8 QVAR(L_LEVELS), WRATIO(L_LEVELS)
75      real*8 KCOEF(4)
76      integer NVAR(L_LEVELS)
77
78  ! temporary variables for multiple aerosol calculation
79  real*8 atemp(L_NLAYRAD,L_NSPECTV)
80  real*8 btemp(L_NLAYRAD,L_NSPECTV)
81  real*8 ctemp(L_NLAYRAD,L_NSPECTV)
82
83!     variables for k in units m^-1
84      real*8 rho, dz
85
86!=======================================================================
87!     Determine the total gas opacity throughout the column, for each
88!     spectral interval, NW, and each Gauss point, NG.
89!     Calculate the continuum opacities, i.e., those that do not depend on
90!     NG, the Gauss index.
91
92taugsurf(:,:) = 0.0
93
94  do K=2,L_LEVELS
95         DPR(K) = PLEV(K)-PLEV(K-1)
96
97
98!         rho  = PLEV(K)/(R*TMID(K))
99         rho  = PMID(K)/(R*TMID(K))
100         dz   = -DPR(k)/(g*rho)
101
102         U(k)   = Cmk*DPR(k)
103
104
105
106
107         call tpindex(PMID(K),TMID(K),QVAR(K),pfgasref,tgasref,WREFVAR, &
108             LCOEF,MT(K),MP(K),NVAR(K),WRATIO(K))
109
110         do LK=1,4
111            LKCOEF(K,LK) = LCOEF(LK)
112         end do
113  end do                    ! levels
114
115  ! Spectral dependance of aerosol absorption
116  do iaer=1,naerkind
117     do NW=1,L_NSPECTV
118        do K=2,L_LEVELS
119           TAEROS(K,NW,IAER) = TAUAERO(K,IAER) * QXVAER(K,NW,IAER)
120        end do
121     end do
122  end do
123
124  ! Rayleigh scattering
125  do NW=1,L_NSPECTV
126     do K=2,L_LEVELS
127        TRAY(K,NW)   = TAURAY(NW) * DPR(K)
128     end do
129  end do
130
131!     we ignore K=1... hope this is ok...
132  do K=2,L_LEVELS
133
134     do NW=1,L_NSPECTV
135
136        TRAYAER = TRAY(K,NW)
137        !     TRAYAER is Tau RAYleigh scattering, plus AERosol opacity
138        do iaer=1,naerkind
139           TRAYAER = TRAYAER + TAEROS(K,NW,IAER)
140        end do
141
142        do NG=1,L_NGAUSS-1
143
144           ! Now compute TAUGAS
145
146           ! Interpolate between water mixing ratios
147           ! WRATIO = 0.0 if the requested water amount is equal to, or outside the
148           ! the water data range
149
150           if (L_REFVAR.eq.1)then ! added by RW for special no variable case
151
152                  KCOEF(1) = GASV(MT(K),MP(K),1,NW,NG)
153                  KCOEF(2) = GASV(MT(K),MP(K)+1,1,NW,NG)
154                  KCOEF(3) = GASV(MT(K)+1,MP(K)+1,1,NW,NG)
155                  KCOEF(4) = GASV(MT(K)+1,MP(K),1,NW,NG)
156
157           else
158
159               KCOEF(1) = GASV(MT(K),MP(K),NVAR(K),NW,NG) + WRATIO(K)*  &
160                   (GASV(MT(K),MP(K),NVAR(K)+1,NW,NG) -                 &
161                   GASV(MT(K),MP(K),NVAR(K),NW,NG))
162
163               KCOEF(2) = GASV(MT(K),MP(K)+1,NVAR(K),NW,NG) + WRATIO(K)* &
164                   (GASV(MT(K),MP(K)+1,NVAR(K)+1,NW,NG) -                &
165                   GASV(MT(K),MP(K)+1,NVAR(K),NW,NG))
166
167             KCOEF(3) = GASV(MT(K)+1,MP(K)+1,NVAR(K),NW,NG) + WRATIO(K)* &
168                   (GASV(MT(K)+1,MP(K)+1,NVAR(K)+1,NW,NG) -  &
169                   GASV(MT(K)+1,MP(K)+1,NVAR(K),NW,NG))
170
171               KCOEF(4) = GASV(MT(K)+1,MP(K),NVAR(K),NW,NG) + WRATIO(K)* &
172                   (GASV(MT(K)+1,MP(K),NVAR(K)+1,NW,NG) -  &
173                   GASV(MT(K)+1,MP(K),NVAR(K),NW,NG))
174           endif
175
176!     Interpolate the gaseous k-coefficients to the requested T,P values
177               
178           ANS = LKCOEF(K,1)*KCOEF(1) + LKCOEF(K,2)*KCOEF(2) + &
179                   LKCOEF(K,3)*KCOEF(3) + LKCOEF(K,4)*KCOEF(4)
180
181           TAUGAS          = U(k)*ANS
182
183           TAUGSURF(NW,NG) = TAUGSURF(NW,NG) + TAUGAS
184           DTAUKV(K,nw,ng) = TAUGAS + TRAYAER
185!               write(21,*) 'TB17 taugas',K,NW,ng,TAUGAS
186
187
188        end do
189
190
191!     Now fill in the "clear" part of the spectrum (NG = L_NGAUSS),
192!     which holds continuum opacity only
193
194        NG = L_NGAUSS
195        DTAUKV(K,nw,ng) = TRAYAER ! Scattering
196
197     end do
198  end do
199
200!=======================================================================
201!     Now the full treatment for the layers, where besides the opacity
202!     we need to calculate the scattering albedo and asymmetry factors
203      !TAUAEROLK(:,:,:) = 1.e-20 ! TB17
204
205  do iaer=1,naerkind
206    DO NW=1,L_NSPECTV
207      DO K=2,L_LEVELS
208           TAUAEROLK(K,NW,IAER) = TAUAERO(K,IAER) * QSVAER(K,NW,IAER) ! effect of scattering albedo
209           !TAUAEROLK(K,NW,IAER) = max(TAUAEROLK(K,NW,IAER),1.e-20) ! TB17
210      end do
211    ENDDO
212  ENDDO
213  !print*, 'TBbug  TAUAEROLK =', TAUAEROLK
214
215  DO NW=1,L_NSPECTV
216     DO L=1,L_NLAYRAD-1
217        K              = 2*L+1
218        atemp(L,NW) = SUM(GVAER(K,NW,1:naerkind) * TAUAEROLK(K,NW,1:naerkind))+SUM(GVAER(K+1,NW,1:naerkind) * TAUAEROLK(K+1,NW,1:naerkind))
219        btemp(L,NW) = SUM(TAUAEROLK(K,NW,1:naerkind)) + SUM(TAUAEROLK(K+1,NW,1:naerkind))
220        ctemp(L,NW) = btemp(L,NW) + 0.9999*(TRAY(K,NW) + TRAY(K+1,NW))  ! JVO 2017 : does this 0.999 is really meaningful ?
221        btemp(L,NW) = btemp(L,NW) + TRAY(K,NW) + TRAY(K+1,NW)
222        COSBV(L,NW,1:L_NGAUSS) = atemp(L,NW)/btemp(L,NW)
223     END DO ! L vertical loop
224     
225     ! Last level
226     L           = L_NLAYRAD
227     K           = 2*L+1
228     atemp(L,NW) = SUM(GVAER(K,NW,1:naerkind) * TAUAEROLK(K,NW,1:naerkind))
229     btemp(L,NW) = SUM(TAUAEROLK(K,NW,1:naerkind))
230     ctemp(L,NW) = btemp(L,NW) + 0.9999*TRAY(K,NW) ! JVO 2017 : does this 0.999 is really meaningful ?
231     btemp(L,NW) = btemp(L,NW) + TRAY(K,NW)
232     COSBV(L,NW,1:L_NGAUSS) = atemp(L,NW)/btemp(L,NW)
233     
234     
235  END DO                    ! NW spectral loop
236
237  DO NG=1,L_NGAUSS
238    DO NW=1,L_NSPECTV
239     DO L=1,L_NLAYRAD-1
240
241        K              = 2*L+1
242        DTAUV(L,nw,ng) = DTAUKV(K,NW,NG) + DTAUKV(K+1,NW,NG)
243        WBARV(L,nw,ng) = ctemp(L,NW) / DTAUV(L,nw,ng)
244
245      END DO ! L vertical loop
246
247        ! Last level
248
249        L              = L_NLAYRAD
250        K              = 2*L+1
251        DTAUV(L,nw,ng) = DTAUKV(K,NW,NG)
252
253        WBARV(L,NW,NG) = ctemp(L,NW) / DTAUV(L,NW,NG)
254        !print*, 'TB22 : WBARV(L)=',WBARV(L,NW,NG),NW,NG
255        !print*, 'TB22 : ctemp(L)=',ctemp(L,NW),NW
256        !print*, 'TB22 : dtauv(L)=',DTAUV(L,NW,NG),NW,NG
257     END DO                 ! NW spectral loop
258  END DO                    ! NG Gauss loop
259
260  ! Total extinction optical depths
261
262  DO NG=1,L_NGAUSS       ! full gauss loop
263     DO NW=1,L_NSPECTV       
264        TAUV(1,NW,NG)=0.0D0
265        DO L=1,L_NLAYRAD
266           TAUV(L+1,NW,NG)=TAUV(L,NW,NG)+DTAUV(L,NW,NG)
267        END DO
268
269        TAUCUMV(1,NW,NG)=0.0D0
270        DO K=2,L_LEVELS
271           TAUCUMV(K,NW,NG)=TAUCUMV(K-1,NW,NG)+DTAUKV(K,NW,NG)
272        END DO
273     END DO           
274  END DO                 ! end full gauss loop
275
276
277  return
278
279
280end subroutine optcv
Note: See TracBrowser for help on using the repository browser.