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

Last change on this file since 1526 was 1397, checked in by milmd, 10 years ago

In LMDZ.GENERIC replacement of all phystd .h files by module files.

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