source: trunk/LMDZ.GENERIC/libf/phystd/optci.F90 @ 880

Last change on this file since 880 was 878, checked in by jleconte, 12 years ago

12/02/2013 == JL

  • Follows previous commit by Aymeric about bilinear interpolations:
    • Extended to all existing continua
    • generalized bilinearbig to work for various size inputs
  • because N2 and H2O continua databases are smaller, improvement around 15% for

an earth case.

File size: 13.3 KB
Line 
1subroutine optci(PLEV,TLEV,DTAUI,TAUCUMI,      &
2     QXIAER,QSIAER,GIAER,COSBI,WBARI,TAUAERO,  &
3     TMID,PMID,TAUGSURF,QVAR,MUVAR)
4
5  use radinc_h
6  use radcommon_h, only: gasi,tlimit,wrefVAR,Cmk,tgasref,pfgasref,wnoi,scalep,indi
7  use gases_h
8  implicit none
9
10  !==================================================================
11  !     
12  !     Purpose
13  !     -------
14  !     Calculates longwave optical constants at each level. For each
15  !     layer and spectral interval in the IR it calculates WBAR, DTAU
16  !     and COSBAR. For each level it calculates TAU.
17  !     
18  !     TAUI(L,LW) is the cumulative optical depth at level L (or alternatively
19  !     at the *bottom* of layer L), LW is the spectral wavelength interval.
20  !     
21  !     TLEV(L) - Temperature at the layer boundary (i.e., level)
22  !     PLEV(L) - Pressure at the layer boundary (i.e., level)
23  !
24  !     Authors
25  !     -------
26  !     Adapted from the NASA Ames code by R. Wordsworth (2009)
27  !     
28  !==================================================================
29
30#include "comcstfi.h"
31#include "callkeys.h"
32
33
34  real*8 DTAUI(L_NLAYRAD,L_NSPECTI,L_NGAUSS)
35  real*8 DTAUKI(L_LEVELS+1,L_NSPECTI,L_NGAUSS)
36  real*8 TAUI(L_NLEVRAD,L_NSPECTI,L_NGAUSS)
37  real*8 TAUCUMI(L_LEVELS,L_NSPECTI,L_NGAUSS)
38  real*8 PLEV(L_LEVELS)
39  real*8 TLEV(L_LEVELS)
40  real*8 TMID(L_LEVELS), PMID(L_LEVELS)
41  real*8 COSBI(L_NLAYRAD,L_NSPECTI,L_NGAUSS)
42  real*8 WBARI(L_NLAYRAD,L_NSPECTI,L_NGAUSS)
43
44  ! for aerosols
45  real*8  QXIAER(L_LEVELS+1,L_NSPECTI,NAERKIND)
46  real*8  QSIAER(L_LEVELS+1,L_NSPECTI,NAERKIND)
47  real*8  GIAER(L_LEVELS+1,L_NSPECTI,NAERKIND)
48  real*8  TAUAERO(L_LEVELS+1,NAERKIND)
49  real*8  TAUAEROLK(L_LEVELS+1,L_NSPECTI,NAERKIND)
50  real*8  TAEROS(L_LEVELS,L_NSPECTI,NAERKIND)
51
52  integer L, NW, NG, K, LK, IAER
53  integer MT(L_LEVELS), MP(L_LEVELS), NP(L_LEVELS)
54  real*8  ANS, TAUGAS
55  real*8  DPR(L_LEVELS), U(L_LEVELS)
56  real*8  LCOEF(4), LKCOEF(L_LEVELS,4)
57
58  real*8 taugsurf(L_NSPECTI,L_NGAUSS-1)
59  real*8 DCONT
60  double precision wn_cont, p_cont, p_air, T_cont, dtemp, dtempc
61  double precision p_cross
62
63  ! variable species mixing ratio variables
64  real*8  QVAR(L_LEVELS), WRATIO(L_LEVELS), MUVAR(L_LEVELS)
65  real*8  KCOEF(4)
66  integer NVAR(L_LEVELS)
67
68  ! temporary variables for multiple aerosol calculation
69  real*8 atemp, btemp
70
71  ! variables for k in units m^-1
72  real*8 dz(L_LEVELS)
73  !real*8 rho !! see test below
74
75  integer igas, jgas
76
77  integer interm
78
79  !--- Kasting's CIA ----------------------------------------
80  !real*8, parameter :: Ci(L_NSPECTI)=[                         &
81  !     3.8E-5, 1.2E-5, 2.8E-6, 7.6E-7, 4.5E-7, 2.3E-7,    &
82  !     5.4E-7, 1.6E-6, 0.0,                               &
83  !     0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0,            &
84  !     0.0, 4.0E-7, 4.0E-6, 1.4E-5,    &
85  !     1.0E-5, 1.2E-6, 2.0E-7, 5.0E-8, 3.0E-8, 0.0 ]
86  !real*8, parameter :: Ti(L_NSPECTI)=[ -2.2, -1.9,             &
87  !     -1.7, -1.7, -1.7, -1.7, -1.7, -1.7,                &
88  !     0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, &
89  !     -1.7,-1.7,-1.7,-1.7,-1.7,-1.7,-1.7, -1.7,0.0 ]
90  !----------------------------------------------------------
91
92  !! AS: to save time in computing continuum (see bilinearbig)
93  IF (.not.ALLOCATED(indi)) THEN
94      ALLOCATE(indi(L_NSPECTI,ngasmx,ngasmx))
95      indi = -9999  ! this initial value means "to be calculated"
96  ENDIF
97
98  !=======================================================================
99  !     Determine the total gas opacity throughout the column, for each
100  !     spectral interval, NW, and each Gauss point, NG.
101
102  taugsurf(:,:) = 0.0
103  dpr(:)        = 0.0
104  lkcoef(:,:)   = 0.0
105
106  do K=2,L_LEVELS
107     DPR(k) = PLEV(K)-PLEV(K-1)
108
109     !--- Kasting's CIA ----------------------------------------
110     !dz(k)=dpr(k)*189.02*TMID(K)/(0.03720*PMID(K))
111     ! this is CO2 path length (in cm) as written by Francois
112     ! delta_z = delta_p * R_specific * T / (g * P)
113     ! But Kasting states that W is in units of _atmosphere_ cm
114     ! So we do
115     !dz(k)=dz(k)*(PMID(K)/1013.25)
116     !dz(k)=dz(k)/100.0 ! in m for SI calc
117     !----------------------------------------------------------
118
119     ! if we have continuum opacities, we need dz
120     if(kastprof)then
121        dz(k) = dpr(k)*(1000.0*8.3145/muvar(k))*TMID(K)/(g*PMID(K))
122        U(k)  = (Cmk*mugaz/(muvar(k)))*DPR(k)
123     else
124        dz(k) = dpr(k)*R*TMID(K)/(g*PMID(K))
125        U(k)  = Cmk*DPR(k)    ! only Cmk line in optci.F
126     endif
127
128     call tpindex(PMID(K),TMID(K),QVAR(K),pfgasref,tgasref,WREFVAR, &
129          LCOEF,MT(K),MP(K),NVAR(K),WRATIO(K))
130
131     do LK=1,4
132        LKCOEF(K,LK) = LCOEF(LK)
133     end do
134
135
136     DO NW=1,L_NSPECTI
137        do iaer=1,naerkind
138           TAEROS(K,NW,IAER) = TAUAERO(K,IAER) * QXIAER(K,NW,IAER)
139        end do
140     END DO
141  end do                    ! levels
142
143  do K=2,L_LEVELS
144
145     do NW=1,L_NSPECTI
146
147        DCONT = 0.0 ! continuum absorption
148
149        if(continuum.and.(.not.graybody))then
150           ! include continua if necessary
151           wn_cont = dble(wnoi(nw))
152           T_cont  = dble(TMID(k))
153           do igas=1,ngasmx
154
155              if(gfrac(igas).eq.-1)then ! variable
156                 p_cont  = dble(PMID(k)*scalep*QVAR(k)) ! qvar = mol/mol
157              else
158                 p_cont  = dble(PMID(k)*scalep*gfrac(igas)*(1.-QVAR(k)))
159              endif
160
161              dtemp=0.0
162              if(igas.eq.igas_N2)then
163
164                 interm = indi(nw,igas,igas)
165                 call interpolateN2N2(wn_cont,T_cont,p_cont,dtemp,.false.,interm)
166                 indi(nw,igas,igas) = interm
167
168              elseif(igas.eq.igas_H2)then
169
170                 ! first do self-induced absorption
171                 interm = indi(nw,igas,igas)
172                 call interpolateH2H2(wn_cont,T_cont,p_cont,dtemp,.false.,interm)
173                 indi(nw,igas,igas) = interm
174
175                 ! then cross-interactions with other gases
176                 do jgas=1,ngasmx
177                    p_cross = dble(PMID(k)*scalep*gfrac(jgas)*(1.-QVAR(k)))
178                    dtempc  = 0.0
179                    if(jgas.eq.igas_N2)then
180                       interm = indi(nw,igas,jgas)
181                       call interpolateN2H2(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm)
182                       indi(nw,igas,jgas) = interm
183                    elseif(jgas.eq.igas_He)then
184                       interm = indi(nw,igas,jgas)
185                       call interpolateH2He(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm)
186                       indi(nw,igas,jgas) = interm
187                    endif
188                    dtemp = dtemp + dtempc
189                 enddo
190
191              elseif(igas.eq.igas_H2O.and.T_cont.gt.200.0)then
192
193                 p_air = dble(PMID(k)*scalep) - p_cont ! note assumes background is air!
194                 if(H2Ocont_simple)then
195                    call interpolateH2Ocont_PPC(wn_cont,T_cont,p_cont,p_air,dtemp,.false.)
196                 else
197                    interm = indi(nw,igas,igas)
198                    call interpolateH2Ocont_CKD(wn_cont,T_cont,p_cont,p_air,dtemp,.false.,interm)
199                    indi(nw,igas,igas) = interm
200                 endif
201
202              endif
203
204              DCONT = DCONT + dtemp
205
206           enddo
207
208           ! Oobleck test
209           !rho = PMID(k)*scalep / (TMID(k)*286.99)
210           !if(WNOI(nw).gt.300.0 .and. WNOI(nw).lt.500.0)then
211           !   DCONT = rho * 0.125 * 4.6e-4
212           !elseif(WNOI(nw).gt.500.0 .and. WNOI(nw).lt.700.0)then
213           !   DCONT = 1000*dpr(k) * 1.0 * 4.6e-4 / g
214           !   DCONT = rho * 1.0 * 4.6e-4
215           !elseif(WNOI(nw).gt.700.0 .and. WNOI(nw).lt.900.0)then
216           !   DCONT = rho * 0.125 * 4.6e-4
217           !endif
218
219           DCONT = DCONT*dz(k)
220
221        endif
222
223        ! RW 7/3/12: already done above
224        !if(.not.Continuum)then
225        !   DCONT=0.0
226        !endif
227
228        !--- Kasting's CIA ----------------------------------------
229        !DCO2   = dz(k)*Ci(nw)*(1.2859*PMID(k)/1000.0)*(TMID(k)/300.)**Ti(nw)
230        !DCO2 = 130*Ci(nw)*(pmid(k)/1013.25)**2*(tmid(k)/300.)**Ti(nw) * dz(k)
231        ! these two have been verified to give the same results
232        !----------------------------------------------------------
233
234        do ng=1,L_NGAUSS-1
235
236           ! Now compute TAUGAS
237
238           ! Interpolate between water mixing ratios
239           ! WRATIO = 0.0 if the requested water amount is equal to, or outside the
240           ! the water data range
241
242           if(L_REFVAR.eq.1)then ! added by RW for special no variable case
243              KCOEF(1) = GASI(MT(K),MP(K),1,NW,NG)
244              KCOEF(2) = GASI(MT(K),MP(K)+1,1,NW,NG)
245              KCOEF(3) = GASI(MT(K)+1,MP(K)+1,1,NW,NG)
246              KCOEF(4) = GASI(MT(K)+1,MP(K),1,NW,NG)
247           else
248
249              KCOEF(1) = GASI(MT(K),MP(K),NVAR(K),NW,NG) + WRATIO(K)*     &
250                   (GASI(MT(K),MP(K),NVAR(K)+1,NW,NG) -                   &
251                   GASI(MT(K),MP(K),NVAR(K),NW,NG))
252
253              KCOEF(2) = GASI(MT(K),MP(K)+1,NVAR(K),NW,NG) + WRATIO(K)*   &
254                   (GASI(MT(K),MP(K)+1,NVAR(K)+1,NW,NG) -                 &
255                   GASI(MT(K),MP(K)+1,NVAR(K),NW,NG))
256
257              KCOEF(3) = GASI(MT(K)+1,MP(K)+1,NVAR(K),NW,NG) + WRATIO(K)* &
258                   (GASI(MT(K)+1,MP(K)+1,NVAR(K)+1,NW,NG) -               &
259                   GASI(MT(K)+1,MP(K)+1,NVAR(K),NW,NG))
260
261              KCOEF(4) = GASI(MT(K)+1,MP(K),NVAR(K),NW,NG) + WRATIO(K)*   &
262                   (GASI(MT(K)+1,MP(K),NVAR(K)+1,NW,NG) -                 &
263                   GASI(MT(K)+1,MP(K),NVAR(K),NW,NG))
264
265           endif
266
267           ! Interpolate the gaseous k-coefficients to the requested T,P values
268
269           ANS = LKCOEF(K,1)*KCOEF(1) + LKCOEF(K,2)*KCOEF(2) +            &
270                LKCOEF(K,3)*KCOEF(3) + LKCOEF(K,4)*KCOEF(4)
271
272           TAUGAS  = U(k)*ANS
273
274           TAUGSURF(NW,NG) = TAUGSURF(NW,NG) + TAUGAS + DCONT
275           DTAUKI(K,nw,ng) = TAUGAS &
276                             + DCONT ! For parameterized continuum absorption
277
278           do iaer=1,naerkind
279              DTAUKI(K,nw,ng) = DTAUKI(K,nw,ng) + TAEROS(K,NW,IAER)
280           end do ! a bug was here!
281
282        end do
283
284        ! Now fill in the "clear" part of the spectrum (NG = L_NGAUSS),
285        ! which holds continuum opacity only
286
287        NG              = L_NGAUSS
288        DTAUKI(K,nw,ng) = 0.0 + DCONT ! For parameterized continuum absorption
289
290        do iaer=1,naerkind
291           DTAUKI(K,nw,ng) = DTAUKI(K,nw,ng) +  TAEROS(K,NW,IAER)
292        end do ! a bug was here!
293
294     end do
295  end do
296
297
298  !=======================================================================
299  !     Now the full treatment for the layers, where besides the opacity
300  !     we need to calculate the scattering albedo and asymmetry factors
301
302  do iaer=1,naerkind
303  DO NW=1,L_NSPECTI
304     DO K=2,L_LEVELS+1
305           TAUAEROLK(K,NW,IAER) = TAUAERO(K,IAER)*QSIAER(K,NW,IAER)
306     ENDDO
307  ENDDO
308  end do
309
310  DO NW=1,L_NSPECTI
311     NG = L_NGAUSS
312     DO L=1,L_NLAYRAD
313
314        K              = 2*L+1
315        DTAUI(L,nw,ng) = DTAUKI(K,NW,NG) + DTAUKI(K+1,NW,NG)! + 1.e-50
316
317        atemp = 0.
318        btemp = 0.
319        if(DTAUI(L,NW,NG) .GT. 1.0E-9) then
320           do iaer=1,naerkind
321              atemp = atemp +                                     &
322                   GIAER(K,NW,IAER)   * TAUAEROLK(K,NW,IAER) +    &
323                   GIAER(K+1,NW,IAER) * TAUAEROLK(K+1,NW,IAER)
324              btemp = btemp + TAUAEROLK(K,NW,IAER) + TAUAEROLK(K+1,NW,IAER)
325              !     *                    + 1.e-10
326           end do
327           WBARI(L,nw,ng) = btemp  / DTAUI(L,NW,NG)
328        else
329           WBARI(L,nw,ng) = 0.0D0
330           DTAUI(L,NW,NG) = 1.0E-9
331        endif
332
333        if(btemp .GT. 0.0) then
334           cosbi(L,NW,NG) = atemp/btemp
335        else
336           cosbi(L,NW,NG) = 0.0D0
337        end if
338
339     END DO ! L vertical loop
340
341     ! Now the other Gauss points, if needed.
342
343     DO NG=1,L_NGAUSS-1
344        IF(TAUGSURF(NW,NG) .gt. TLIMIT) THEN
345
346           DO L=1,L_NLAYRAD
347              K              = 2*L+1
348              DTAUI(L,nw,ng) = DTAUKI(K,NW,NG)+DTAUKI(K+1,NW,NG)! + 1.e-50
349
350              btemp = 0.
351              if(DTAUI(L,NW,NG) .GT. 1.0E-9) then
352
353                 do iaer=1,naerkind
354                    btemp = btemp + TAUAEROLK(K,NW,IAER) + TAUAEROLK(K+1,NW,IAER)
355                 end do
356                 WBARI(L,nw,ng) = btemp  / DTAUI(L,NW,NG)
357
358              else
359                 WBARI(L,nw,ng) = 0.0D0
360                 DTAUI(L,NW,NG) = 1.0E-9
361              endif
362
363              cosbi(L,NW,NG) = cosbi(L,NW,L_NGAUSS)
364           END DO ! L vertical loop
365        END IF
366
367     END DO                 ! NG Gauss loop
368  END DO                    ! NW spectral loop
369
370  ! Total extinction optical depths
371
372  DO NW=1,L_NSPECTI       
373     DO NG=1,L_NGAUSS       ! full gauss loop
374        TAUI(1,NW,NG)=0.0D0
375        DO L=1,L_NLAYRAD
376           TAUI(L+1,NW,NG)=TAUI(L,NW,NG)+DTAUI(L,NW,NG)
377        END DO
378
379        TAUCUMI(1,NW,NG)=0.0D0
380        DO K=2,L_LEVELS
381           TAUCUMI(K,NW,NG)=TAUCUMI(K-1,NW,NG)+DTAUKI(K,NW,NG)
382        END DO
383     END DO                 ! end full gauss loop
384  END DO
385
386  ! be aware when comparing with textbook results
387  ! (e.g. Pierrehumbert p. 218) that
388  ! taucumi does not take the <cos theta>=0.5 factor into
389  ! account. It is the optical depth for a vertically
390  ! ascending ray with angle theta = 0.
391
392  !open(127,file='taucum.out')
393  !do nw=1,L_NSPECTI
394  !   write(127,*) taucumi(L_LEVELS,nw,L_NGAUSS)
395  !enddo
396  !close(127)
397
398  return
399
400
401end subroutine optci
402
403
404
Note: See TracBrowser for help on using the repository browser.