source: trunk/LMDZ.GENERIC/libf/phystd/optcv.F90 @ 889

Last change on this file since 889 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.

  • Property svn:executable set to *
File size: 11.3 KB
RevLine 
[716]1SUBROUTINE OPTCV(DTAUV,TAUV,TAUCUMV,PLEV,  &
2     QXVAER,QSVAER,GVAER,WBARV,COSBV,       &
3     TAURAY,TAUAERO,TMID,PMID,TAUGSURF,QVAR,MUVAR)
[253]4
[716]5  use radinc_h
[873]6  use radcommon_h, only: gasv, tlimit, wrefVAR, Cmk, tgasref, pfgasref,wnov,scalep,indv
[716]7  use gases_h
[253]8
[716]9  implicit none
[253]10
[716]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  !     THIS SUBROUTINE SETS THE OPTICAL CONSTANTS IN THE VISUAL 
24  !     IT CALCUALTES FOR EACH LAYER, FOR EACH SPECRAL INTERVAL IN THE VISUAL
25  !     LAYER: WBAR, DTAU, COSBAR
26  !     LEVEL: TAU
27  !     
28  !     TAUV(L,NW,NG) is the cumulative optical depth at the top of radiation code
29  !     layer L. NW is spectral wavelength interval, ng the Gauss point index.
30  !     
31  !     TLEV(L) - Temperature at the layer boundary
32  !     PLEV(L) - Pressure at the layer boundary (i.e. level)
33  !     GASV(NT,NPS,NW,NG) - Visible k-coefficients
34  !     
35  !-------------------------------------------------------------------
[253]36
[873]37#include "comcstfi.h"
[253]38#include "callkeys.h"
39
[873]40
[716]41  real*8 DTAUV(L_NLAYRAD,L_NSPECTV,L_NGAUSS)
42  real*8 DTAUKV(L_LEVELS+1,L_NSPECTV,L_NGAUSS)
43  real*8 TAUV(L_NLEVRAD,L_NSPECTV,L_NGAUSS)
44  real*8 TAUCUMV(L_LEVELS,L_NSPECTV,L_NGAUSS)
45  real*8 PLEV(L_LEVELS)
46  real*8 TMID(L_LEVELS), PMID(L_LEVELS)
47  real*8 COSBV(L_NLAYRAD,L_NSPECTV,L_NGAUSS)
48  real*8 WBARV(L_NLAYRAD,L_NSPECTV,L_NGAUSS)
[253]49
[716]50  ! for aerosols
[873]51  real*8  QXVAER(L_LEVELS+1,L_NSPECTV,NAERKIND)
52  real*8  QSVAER(L_LEVELS+1,L_NSPECTV,NAERKIND)
53  real*8  GVAER(L_LEVELS+1,L_NSPECTV,NAERKIND)
54  real*8  TAUAERO(L_LEVELS+1,NAERKIND)
55  real*8  TAUAEROLK(L_LEVELS+1,L_NSPECTV,NAERKIND)
56  real*8  TAEROS(L_LEVELS,L_NSPECTV,NAERKIND)
[253]57
[873]58  integer L, NW, NG, K, LK, IAER
[716]59  integer MT(L_LEVELS), MP(L_LEVELS), NP(L_LEVELS)
60  real*8  ANS, TAUGAS
[873]61  real*8  TAURAY(L_NSPECTV)
[716]62  real*8  TRAY(L_LEVELS,L_NSPECTV)
[873]63  real*8  TRAYAER
[716]64  real*8  DPR(L_LEVELS), U(L_LEVELS)
65  real*8  LCOEF(4), LKCOEF(L_LEVELS,4)
[253]66
[873]67  real*8 taugsurf(L_NSPECTV,L_NGAUSS-1)
68  real*8 DCONT
69  double precision wn_cont, p_cont, p_air, T_cont, dtemp, dtempc
70  double precision p_cross
[253]71
[716]72  ! variable species mixing ratio variables
[873]73  real*8  QVAR(L_LEVELS), WRATIO(L_LEVELS), MUVAR(L_LEVELS)
74  real*8  KCOEF(4)
[716]75  integer NVAR(L_LEVELS)
[253]76
[716]77  ! temporary variables for multiple aerosol calculation
78  real*8 atemp, btemp, ctemp
[253]79
[716]80  ! variables for k in units m^-1
[873]81  real*8 dz(L_LEVELS)
[253]82
[716]83  integer igas, jgas
[253]84
[873]85  integer interm
86
87  !! AS: to save time in computing continuum (see bilinearbig)
88  IF (.not.ALLOCATED(indv)) THEN
[878]89      ALLOCATE(indv(L_NSPECTV,ngasmx,ngasmx))
[873]90      indv = -9999 ! this initial value means "to be calculated"
91  ENDIF
92
[716]93  !=======================================================================
94  !     Determine the total gas opacity throughout the column, for each
95  !     spectral interval, NW, and each Gauss point, NG.
96  !     Calculate the continuum opacities, i.e., those that do not depend on
97  !     NG, the Gauss index.
[253]98
[716]99  taugsurf(:,:) = 0.0
100  dpr(:)        = 0.0
101  lkcoef(:,:)   = 0.0
[253]102
[716]103  do K=2,L_LEVELS
104     DPR(k) = PLEV(K)-PLEV(K-1)
[253]105
[716]106     ! if we have continuum opacities, we need dz
107     if(kastprof)then
[873]108        dz(k) = dpr(k)*(1000.0*8.3145/muvar(k))*TMID(K)/(g*PMID(K))
[716]109        U(k)  = (Cmk*mugaz/(muvar(k)))*DPR(k)
110     else
111        dz(k) = dpr(k)*R*TMID(K)/(g*PMID(K))
112        U(k)  = Cmk*DPR(k)    ! only Cmk line in optci.F
113     endif
[253]114
[716]115     call tpindex(PMID(K),TMID(K),QVAR(K),pfgasref,tgasref,WREFVAR, &
116          LCOEF,MT(K),MP(K),NVAR(K),WRATIO(K))
[253]117
[716]118     do LK=1,4
119        LKCOEF(K,LK) = LCOEF(LK)
120     end do
[253]121
[873]122
[716]123     DO NW=1,L_NSPECTV
124        do iaer=1,naerkind
[873]125           TAEROS(K,NW,IAER) = TAUAERO(K,IAER) * QXVAER(K,NW,IAER)
[716]126        end do
[873]127        TRAY(K,NW)   = TAURAY(NW) * DPR(K)
[716]128     END DO
[873]129  end do                    ! levels
[253]130
[716]131  !     we ignore K=1...
132  do K=2,L_LEVELS
[873]133
[716]134     do NW=1,L_NSPECTV
[253]135
[716]136        TRAYAER = TRAY(K,NW)
[873]137        !     TRAYAER is Tau RAYleigh scattering, plus AERosol opacity
[716]138        do iaer=1,naerkind
139           TRAYAER = TRAYAER + TAEROS(K,NW,IAER)
140        end do
[253]141
[716]142        DCONT = 0.0 ! continuum absorption
[253]143
[873]144        if(continuum.and.(.not.graybody).and.callgasvis)then
[716]145           ! include continua if necessary
146           wn_cont = dble(wnov(nw))
147           T_cont  = dble(TMID(k))
148           do igas=1,ngasmx
[305]149
[716]150              if(gfrac(igas).eq.-1)then ! variable
151                 p_cont  = dble(PMID(k)*scalep*QVAR(k)) ! qvar = mol/mol
152              else
153                 p_cont  = dble(PMID(k)*scalep*gfrac(igas)*(1.-QVAR(k)))
154              endif
[305]155
[716]156              dtemp=0.0
157              if(igas.eq.igas_N2)then
[253]158
[878]159                 interm = indv(nw,igas,igas)
160!                 call interpolateN2N2(wn_cont,T_cont,p_cont,dtemp,.false.,interm)
161                 indv(nw,igas,igas) = interm
[716]162                 ! only goes to 500 cm^-1, so unless we're around a cold brown dwarf, this is irrelevant in the visible
[253]163
[716]164              elseif(igas.eq.igas_H2)then
[253]165
[716]166                 ! first do self-induced absorption
[878]167                 interm = indv(nw,igas,igas)
[873]168                 call interpolateH2H2(wn_cont,T_cont,p_cont,dtemp,.false.,interm)
[878]169                 indv(nw,igas,igas) = interm
[253]170
[716]171                 ! then cross-interactions with other gases
172                 do jgas=1,ngasmx
173                    p_cross = dble(PMID(k)*scalep*gfrac(jgas)*(1.-QVAR(k)))
[873]174                    dtempc  = 0.0
175                    if(jgas.eq.igas_N2)then
[878]176                       interm = indv(nw,igas,jgas)
177                       call interpolateN2H2(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm)
178                       indv(nw,igas,jgas) = interm
[716]179                       ! should be irrelevant in the visible
180                    elseif(jgas.eq.igas_He)then
[878]181                       interm = indv(nw,igas,jgas)
[873]182                       call interpolateH2He(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm)
[878]183                       indv(nw,igas,jgas) = interm
[716]184                    endif
[873]185                    dtemp = dtemp + dtempc
[716]186                 enddo
[253]187
[716]188              elseif(igas.eq.igas_H2O.and.T_cont.gt.200.0)then
[253]189
[716]190                 p_air = dble(PMID(k)*scalep) - p_cont ! note assumes background is air!
191                 if(H2Ocont_simple)then
192                    call interpolateH2Ocont_PPC(wn_cont,T_cont,p_cont,p_air,dtemp,.false.)
193                 else
[878]194                    interm = indv(nw,igas,igas)
195                    call interpolateH2Ocont_CKD(wn_cont,T_cont,p_cont,p_air,dtemp,.false.,interm)
196                    indv(nw,igas,igas) = interm
[716]197                 endif
[253]198
[716]199              endif
[253]200
[716]201              DCONT = DCONT + dtemp
[253]202
[716]203           enddo
[253]204
[873]205           DCONT = DCONT*dz(k)
206
[716]207        endif
[253]208
[873]209        do ng=1,L_NGAUSS-1
[305]210
[873]211           ! Now compute TAUGAS
[253]212
[873]213           ! Interpolate between water mixing ratios
214           ! WRATIO = 0.0 if the requested water amount is equal to, or outside the
215           ! the water data range
216
217           if(L_REFVAR.eq.1)then ! added by RW for special no variable case
[716]218              KCOEF(1) = GASV(MT(K),MP(K),1,NW,NG)
219              KCOEF(2) = GASV(MT(K),MP(K)+1,1,NW,NG)
220              KCOEF(3) = GASV(MT(K)+1,MP(K)+1,1,NW,NG)
221              KCOEF(4) = GASV(MT(K)+1,MP(K),1,NW,NG)
222           else
[873]223
[716]224              KCOEF(1) = GASV(MT(K),MP(K),NVAR(K),NW,NG) + WRATIO(K)*    &
225                   (GASV(MT(K),MP(K),NVAR(K)+1,NW,NG) -                  &
226                   GASV(MT(K),MP(K),NVAR(K),NW,NG))
[253]227
[716]228              KCOEF(2) = GASV(MT(K),MP(K)+1,NVAR(K),NW,NG) + WRATIO(K)*  &
229                   (GASV(MT(K),MP(K)+1,NVAR(K)+1,NW,NG) -                &
230                   GASV(MT(K),MP(K)+1,NVAR(K),NW,NG))
[253]231
[716]232              KCOEF(3) = GASV(MT(K)+1,MP(K)+1,NVAR(K),NW,NG) + WRATIO(K)*&
233                   (GASV(MT(K)+1,MP(K)+1,NVAR(K)+1,NW,NG) -              &
234                   GASV(MT(K)+1,MP(K)+1,NVAR(K),NW,NG))
[253]235
[716]236              KCOEF(4) = GASV(MT(K)+1,MP(K),NVAR(K),NW,NG) + WRATIO(K)*  &
237                   (GASV(MT(K)+1,MP(K),NVAR(K)+1,NW,NG) -                &
238                   GASV(MT(K)+1,MP(K),NVAR(K),NW,NG))
[873]239
[716]240           endif
[253]241
[873]242           ! Interpolate the gaseous k-coefficients to the requested T,P values
[253]243
[873]244           ANS = LKCOEF(K,1)*KCOEF(1) + LKCOEF(K,2)*KCOEF(2) +            &
[716]245                LKCOEF(K,3)*KCOEF(3) + LKCOEF(K,4)*KCOEF(4)
[253]246
[873]247           TAUGAS  = U(k)*ANS
[253]248
[716]249           TAUGSURF(NW,NG) = TAUGSURF(NW,NG) + TAUGAS + DCONT
[873]250           DTAUKV(K,nw,ng) = TAUGAS &
251                             + TRAYAER & ! TRAYAER includes all scattering contributions
252                             + DCONT ! For parameterized continuum aborption
[253]253
[716]254        end do
[253]255
[873]256        ! Now fill in the "clear" part of the spectrum (NG = L_NGAUSS),
257        ! which holds continuum opacity only
[253]258
[873]259        NG              = L_NGAUSS
260        DTAUKV(K,nw,ng) = TRAY(K,NW) + DCONT ! For parameterized continuum absorption
[253]261
[716]262        do iaer=1,naerkind
[873]263           DTAUKV(K,nw,ng) = DTAUKV(K,nw,ng) +  TAEROS(K,NW,IAER)
[716]264        end do ! a bug was here!
[253]265
[716]266     end do
267  end do
[253]268
269
[716]270  !=======================================================================
271  !     Now the full treatment for the layers, where besides the opacity
272  !     we need to calculate the scattering albedo and asymmetry factors
[253]273
[873]274  do iaer=1,naerkind
[716]275  DO NW=1,L_NSPECTV
[873]276     DO K=2,L_LEVELS   ! AS: shouldn't this be L_LEVELS+1 ? (see optci)
[716]277           TAUAEROLK(K,NW,IAER) = TAUAERO(K,IAER) * QSVAER(K,NW,IAER)
278     ENDDO
279  ENDDO
[873]280  end do
[253]281
[716]282  DO NW=1,L_NSPECTV
[873]283    DO NG=1,L_NGAUSS
284     DO L=1,L_NLAYRAD-1
[253]285
[873]286        K              = 2*L+1
287        DTAUV(L,nw,ng) = DTAUKV(K,NW,NG) + DTAUKV(K+1,NW,NG)
[253]288
[873]289        atemp = 0.
290        btemp = TRAY(K,NW) + TRAY(K+1,NW)
291        ctemp=0.9999*(TRAY(K,NW) + TRAY(K+1,NW))
[716]292           do iaer=1,naerkind
293              atemp = atemp +                                     &
294                   GVAER(K,NW,IAER)   * TAUAEROLK(K,NW,IAER) +    &
295                   GVAER(K+1,NW,IAER) * TAUAEROLK(K+1,NW,IAER)
[873]296              btemp = btemp + TAUAEROLK(K,NW,IAER) + TAUAEROLK(K+1,NW,IAER)
297              ctemp = ctemp + TAUAEROLK(K,NW,IAER) + TAUAEROLK(K+1,NW,IAER)
[716]298           end do
[873]299           WBARV(L,nw,ng) = ctemp / DTAUV(L,nw,ng)
[716]300           COSBV(L,NW,NG) = atemp/btemp
[253]301
[873]302      END DO ! L vertical loop
[253]303
[716]304        !     No vertical averaging on bottom layer
[253]305
[716]306        L              = L_NLAYRAD
307        K              = 2*L+1
308        DTAUV(L,nw,ng) = DTAUKV(K,NW,NG)
309
310        atemp=0.
311        btemp=TRAY(K,NW)
312        ctemp=0.9999*TRAY(K,NW)
313        do iaer=1,naerkind
314           atemp = atemp + GVAER(K,NW,IAER) * TAUAEROLK(K,NW,IAER)
315           btemp = btemp + TAUAEROLK(K,NW,IAER)
316           ctemp = ctemp + TAUAEROLK(K,NW,IAER)
317        end do
318        COSBV(L,NW,NG) = atemp/btemp
319        WBARV(L,nw,ng) = ctemp/DTAUV(L,nw,ng)
320
[873]321     END DO                 ! NG Gauss loop
[716]322  END DO                    ! NW spectral loop
323
324  ! Total extinction optical depths
325
[873]326  DO NW=1,L_NSPECTV       
[716]327     DO NG=1,L_NGAUSS       ! full gauss loop
328        TAUV(1,NW,NG)=0.0D0
329        DO L=1,L_NLAYRAD
330           TAUV(L+1,NW,NG)=TAUV(L,NW,NG)+DTAUV(L,NW,NG)
331        END DO
332
333        TAUCUMV(1,NW,NG)=0.0D0
334        DO K=2,L_LEVELS
335           TAUCUMV(K,NW,NG)=TAUCUMV(K-1,NW,NG)+DTAUKV(K,NW,NG)
336        END DO
337     END DO                 ! end full gauss loop
338  END DO
339
340
[873]341  return
342
343
344end subroutine optcv
Note: See TracBrowser for help on using the repository browser.