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

Last change on this file since 1990 was 1987, checked in by jleconte, 7 years ago

28/08/2018 == JL

Start a series of commits to change the upper boundary conditions in the radiative transfer to solve some issues with the last two layers.
It seems to be good to have aerosols in the first "radiative layer" of the gcm in the IR but visible does not handle very well diffusion in first layer.
Tauaero and tauray are set to 0 (a small value for rayleigh because the code crashes otherwise) in the 4 first semilayers in optcv, but not optci.
This solves random variations of the sw heating at the model top.

  • Property svn:executable set to *
File size: 12.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 
[1715]26  !     IT CALCULATES FOR EACH LAYER, FOR EACH SPECTRAL INTERVAL IN THE VISUAL
[716]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)
[1715]41  real*8 DTAUKV(L_LEVELS,L_NSPECTV,L_NGAUSS)
[716]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
[1715]50  real*8  QXVAER(L_LEVELS,L_NSPECTV,NAERKIND)
51  real*8  QSVAER(L_LEVELS,L_NSPECTV,NAERKIND)
52  real*8  GVAER(L_LEVELS,L_NSPECTV,NAERKIND)
53  real*8  TAUAERO(L_LEVELS,NAERKIND)
54  real*8  TAUAEROLK(L_LEVELS,L_NSPECTV,NAERKIND)
[873]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)
62  real*8  DPR(L_LEVELS), U(L_LEVELS)
63  real*8  LCOEF(4), LKCOEF(L_LEVELS,4)
[253]64
[873]65  real*8 taugsurf(L_NSPECTV,L_NGAUSS-1)
[918]66  real*8 DCONT,DAERO
[1715]67  real*8 DRAYAER
[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)
[1725]75 
76  ! temporary variables to reduce memory access time to gasv
77  real*8 tmpk(2,2)
78  real*8 tmpkvar(2,2,2)
[253]79
[716]80  ! temporary variables for multiple aerosol calculation
[918]81  real*8 atemp(L_NLAYRAD,L_NSPECTV)
82  real*8 btemp(L_NLAYRAD,L_NSPECTV)
83  real*8 ctemp(L_NLAYRAD,L_NSPECTV)
[253]84
[716]85  ! variables for k in units m^-1
[873]86  real*8 dz(L_LEVELS)
[253]87
[716]88  integer igas, jgas
[253]89
[873]90  integer interm
91
92  !! AS: to save time in computing continuum (see bilinearbig)
93  IF (.not.ALLOCATED(indv)) THEN
[878]94      ALLOCATE(indv(L_NSPECTV,ngasmx,ngasmx))
[873]95      indv = -9999 ! this initial value means "to be calculated"
96  ENDIF
97
[716]98  !=======================================================================
99  !     Determine the total gas opacity throughout the column, for each
100  !     spectral interval, NW, and each Gauss point, NG.
101  !     Calculate the continuum opacities, i.e., those that do not depend on
102  !     NG, the Gauss index.
[253]103
[716]104  taugsurf(:,:) = 0.0
105  dpr(:)        = 0.0
106  lkcoef(:,:)   = 0.0
[253]107
[716]108  do K=2,L_LEVELS
109     DPR(k) = PLEV(K)-PLEV(K-1)
[253]110
[716]111     ! if we have continuum opacities, we need dz
112     if(kastprof)then
[1016]113        dz(k) = dpr(k)*(1000.0d0*8.3145d0/muvar(k))*TMID(K)/(g*PMID(K))
114        U(k)  = Cmk*DPR(k)*mugaz/muvar(k)
[716]115     else
[1194]116        dz(k) = dpr(k)*R*TMID(K)/(glat_ig*PMID(K))*mugaz/muvar(k)
[1016]117        U(k)  = Cmk*DPR(k)*mugaz/muvar(k)     ! only Cmk line in optci.F 
118            !JL13 the mugaz/muvar factor takes into account water meanmolecular weight if water is present
[716]119     endif
[253]120
[716]121     call tpindex(PMID(K),TMID(K),QVAR(K),pfgasref,tgasref,WREFVAR, &
122          LCOEF,MT(K),MP(K),NVAR(K),WRATIO(K))
[253]123
[716]124     do LK=1,4
125        LKCOEF(K,LK) = LCOEF(LK)
126     end do
[918]127  end do                    ! levels
[253]128
[1715]129  ! Spectral dependance of aerosol absorption
[1987]130            !JL18 It seems to be good to have aerosols in the first "radiative layer" of the gcm in the IR
131            !   but visible does not handle very well diffusion in first layer.
132            !   The tauaero and tauray are thus set to 0 (a small value for rayleigh because the code crashes otherwise)
133            !   in the 4 first semilayers in optcv, but not optci.
134            !   This solves random variations of the sw heating at the model top.
[918]135  do iaer=1,naerkind
136     do NW=1,L_NSPECTV
[1987]137        TAEROS(1:4,NW,IAER)=0.d0
138        do K=5,L_LEVELS
[873]139           TAEROS(K,NW,IAER) = TAUAERO(K,IAER) * QXVAER(K,NW,IAER)
[918]140        end do                    ! levels
141     end do
142  end do
[1715]143 
144  ! Rayleigh scattering
[918]145  do NW=1,L_NSPECTV
[1987]146     TRAY(1:4,NW)   = 1d-30
147     do K=5,L_LEVELS
[873]148        TRAY(K,NW)   = TAURAY(NW) * DPR(K)
[918]149     end do                    ! levels
150  end do
151 
[716]152  !     we ignore K=1...
153  do K=2,L_LEVELS
[873]154
[716]155     do NW=1,L_NSPECTV
[253]156
[1715]157        DRAYAER = TRAY(K,NW)
158        !     DRAYAER is Tau RAYleigh scattering, plus AERosol opacity
[716]159        do iaer=1,naerkind
[1715]160           DRAYAER = DRAYAER + TAEROS(K,NW,IAER)
[716]161        end do
[253]162
[716]163        DCONT = 0.0 ! continuum absorption
[253]164
[873]165        if(continuum.and.(.not.graybody).and.callgasvis)then
[716]166           ! include continua if necessary
167           wn_cont = dble(wnov(nw))
168           T_cont  = dble(TMID(k))
169           do igas=1,ngasmx
[305]170
[716]171              if(gfrac(igas).eq.-1)then ! variable
172                 p_cont  = dble(PMID(k)*scalep*QVAR(k)) ! qvar = mol/mol
173              else
174                 p_cont  = dble(PMID(k)*scalep*gfrac(igas)*(1.-QVAR(k)))
175              endif
[305]176
[716]177              dtemp=0.0
178              if(igas.eq.igas_N2)then
[253]179
[878]180                 interm = indv(nw,igas,igas)
181!                 call interpolateN2N2(wn_cont,T_cont,p_cont,dtemp,.false.,interm)
182                 indv(nw,igas,igas) = interm
[716]183                 ! only goes to 500 cm^-1, so unless we're around a cold brown dwarf, this is irrelevant in the visible
[253]184
[716]185              elseif(igas.eq.igas_H2)then
[253]186
[716]187                 ! first do self-induced absorption
[878]188                 interm = indv(nw,igas,igas)
[873]189                 call interpolateH2H2(wn_cont,T_cont,p_cont,dtemp,.false.,interm)
[878]190                 indv(nw,igas,igas) = interm
[253]191
[716]192                 ! then cross-interactions with other gases
193                 do jgas=1,ngasmx
194                    p_cross = dble(PMID(k)*scalep*gfrac(jgas)*(1.-QVAR(k)))
[873]195                    dtempc  = 0.0
196                    if(jgas.eq.igas_N2)then
[878]197                       interm = indv(nw,igas,jgas)
198                       call interpolateN2H2(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm)
199                       indv(nw,igas,jgas) = interm
[716]200                       ! should be irrelevant in the visible
201                    elseif(jgas.eq.igas_He)then
[878]202                       interm = indv(nw,igas,jgas)
[873]203                       call interpolateH2He(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm)
[878]204                       indv(nw,igas,jgas) = interm
[716]205                    endif
[873]206                    dtemp = dtemp + dtempc
[716]207                 enddo
[253]208
[716]209              elseif(igas.eq.igas_H2O.and.T_cont.gt.200.0)then
[253]210
[716]211                 p_air = dble(PMID(k)*scalep) - p_cont ! note assumes background is air!
212                 if(H2Ocont_simple)then
213                    call interpolateH2Ocont_PPC(wn_cont,T_cont,p_cont,p_air,dtemp,.false.)
214                 else
[878]215                    interm = indv(nw,igas,igas)
216                    call interpolateH2Ocont_CKD(wn_cont,T_cont,p_cont,p_air,dtemp,.false.,interm)
217                    indv(nw,igas,igas) = interm
[716]218                 endif
[253]219
[716]220              endif
[253]221
[716]222              DCONT = DCONT + dtemp
[253]223
[716]224           enddo
[253]225
[873]226           DCONT = DCONT*dz(k)
227
[716]228        endif
[253]229
[873]230        do ng=1,L_NGAUSS-1
[305]231
[873]232           ! Now compute TAUGAS
[253]233
[873]234           ! Interpolate between water mixing ratios
235           ! WRATIO = 0.0 if the requested water amount is equal to, or outside the
236           ! the water data range
237
238           if(L_REFVAR.eq.1)then ! added by RW for special no variable case
[1725]239           
240              ! JVO 2017 : added tmpk because the repeated calls to gasi/v increased dramatically
241              ! the execution time of optci/v -> ~ factor 2 on the whole radiative
242              ! transfer on the tested simulations !
243
244              tmpk = GASV(MT(K):MT(K)+1,MP(K):MP(K)+1,1,NW,NG)
245             
246              KCOEF(1) = tmpk(1,1) ! KCOEF(1) = GASV(MT(K),MP(K),1,NW,NG)
247              KCOEF(2) = tmpk(1,2) ! KCOEF(2) = GASV(MT(K),MP(K)+1,1,NW,NG)
248              KCOEF(3) = tmpk(2,2) ! KCOEF(3) = GASV(MT(K)+1,MP(K)+1,1,NW,NG)
249              KCOEF(4) = tmpk(2,1) ! KCOEF(4) = GASV(MT(K)+1,MP(K),1,NW,NG)
250
[716]251           else
[873]252
[1725]253              tmpkvar = GASV(MT(K):MT(K)+1,MP(K):MP(K)+1,NVAR(K):NVAR(K)+1,NW,NG)
[253]254
[1725]255              KCOEF(1) = tmpkvar(1,1,1) + WRATIO(K) *  &
256                        ( tmpkvar(1,1,2)-tmpkvar(1,1,1) )
[253]257
[1725]258              KCOEF(2) = tmpkvar(1,2,1) + WRATIO(K) *  &
259                        ( tmpkvar(1,2,2)-tmpkvar(1,2,1) )
[253]260
[1725]261              KCOEF(3) = tmpkvar(2,2,1) + WRATIO(K) *  &
262                        ( tmpkvar(2,2,2)-tmpkvar(2,2,1) )
263             
264              KCOEF(4) = tmpkvar(2,1,1) + WRATIO(K) *  &
265                        ( tmpkvar(2,1,2)-tmpkvar(2,1,1) )
[873]266
[1725]267
[716]268           endif
[253]269
[873]270           ! Interpolate the gaseous k-coefficients to the requested T,P values
[253]271
[873]272           ANS = LKCOEF(K,1)*KCOEF(1) + LKCOEF(K,2)*KCOEF(2) +            &
[716]273                LKCOEF(K,3)*KCOEF(3) + LKCOEF(K,4)*KCOEF(4)
[253]274
[873]275           TAUGAS  = U(k)*ANS
[253]276
[716]277           TAUGSURF(NW,NG) = TAUGSURF(NW,NG) + TAUGAS + DCONT
[873]278           DTAUKV(K,nw,ng) = TAUGAS &
[1715]279                             + DRAYAER & ! DRAYAER includes all scattering contributions
[873]280                             + DCONT ! For parameterized continuum aborption
[253]281
[716]282        end do
[253]283
[873]284        ! Now fill in the "clear" part of the spectrum (NG = L_NGAUSS),
285        ! which holds continuum opacity only
[253]286
[873]287        NG              = L_NGAUSS
[1715]288        DTAUKV(K,nw,ng) = DRAYAER + DCONT ! Scattering + parameterized continuum absorption
[253]289
[716]290     end do
291  end do
[253]292
293
[716]294  !=======================================================================
295  !     Now the full treatment for the layers, where besides the opacity
296  !     we need to calculate the scattering albedo and asymmetry factors
[253]297
[1987]298            !JL18 It seems to be good to have aerosols in the first "radiative layer" of the gcm in the IR
299            !   but not in the visible
300            !   The tauaero is thus set to 0 in the 4 first semilayers in optcv, but not optci.
301            !   This solves random variations of the sw heating at the model top.
[873]302  do iaer=1,naerkind
[918]303    DO NW=1,L_NSPECTV
[1987]304      TAUAEROLK(1:4,NW,IAER)=0.d0
305      DO K=5,L_LEVELS
[1715]306           TAUAEROLK(K,NW,IAER) = TAUAERO(K,IAER) * QSVAER(K,NW,IAER) ! effect of scattering albedo
[918]307      ENDDO
308    ENDDO
[873]309  end do
[253]310
[716]311  DO NW=1,L_NSPECTV
[919]312     DO L=1,L_NLAYRAD-1
[918]313        K              = 2*L+1
314        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))
315        btemp(L,NW) = SUM(TAUAEROLK(K,NW,1:naerkind)) + SUM(TAUAEROLK(K+1,NW,1:naerkind))
[1715]316        ctemp(L,NW) = btemp(L,NW) + 0.9999*(TRAY(K,NW) + TRAY(K+1,NW))  ! JVO 2017 : does this 0.999 is really meaningful ?
[918]317        btemp(L,NW) = btemp(L,NW) + TRAY(K,NW) + TRAY(K+1,NW)
318        COSBV(L,NW,1:L_NGAUSS) = atemp(L,NW)/btemp(L,NW)
319     END DO ! L vertical loop
[919]320     
[1715]321     ! Last level
322     L           = L_NLAYRAD
323     K           = 2*L+1
324     atemp(L,NW) = SUM(GVAER(K,NW,1:naerkind) * TAUAEROLK(K,NW,1:naerkind))
[919]325     btemp(L,NW) = SUM(TAUAEROLK(K,NW,1:naerkind))
[1715]326     ctemp(L,NW) = btemp(L,NW) + 0.9999*TRAY(K,NW) ! JVO 2017 : does this 0.999 is really meaningful ?
[919]327     btemp(L,NW) = btemp(L,NW) + TRAY(K,NW)
328     COSBV(L,NW,1:L_NGAUSS) = atemp(L,NW)/btemp(L,NW)
329     
330     
[918]331  END DO                    ! NW spectral loop
332
333  DO NG=1,L_NGAUSS
334    DO NW=1,L_NSPECTV
[873]335     DO L=1,L_NLAYRAD-1
[253]336
[873]337        K              = 2*L+1
338        DTAUV(L,nw,ng) = DTAUKV(K,NW,NG) + DTAUKV(K+1,NW,NG)
[918]339        WBARV(L,nw,ng) = ctemp(L,NW) / DTAUV(L,nw,ng)
[253]340
[873]341      END DO ! L vertical loop
[253]342
[1715]343        ! Last level
[253]344
[716]345        L              = L_NLAYRAD
346        K              = 2*L+1
[919]347        DTAUV(L,nw,ng) = DTAUKV(K,NW,NG)
348
349        WBARV(L,NW,NG) = ctemp(L,NW) / DTAUV(L,NW,NG)
[1722]350
[918]351     END DO                 ! NW spectral loop
352  END DO                    ! NG Gauss loop
[716]353
354  ! Total extinction optical depths
355
[918]356  DO NG=1,L_NGAUSS       ! full gauss loop
357     DO NW=1,L_NSPECTV       
[716]358        TAUCUMV(1,NW,NG)=0.0D0
359        DO K=2,L_LEVELS
360           TAUCUMV(K,NW,NG)=TAUCUMV(K-1,NW,NG)+DTAUKV(K,NW,NG)
361        END DO
[1987]362
363        DO L=1,L_NLAYRAD+1
364           TAUV(L,NW,NG)=TAUCUMV(2*L,NW,NG)
365        END DO
[918]366     END DO           
367  END DO                 ! end full gauss loop
[716]368
369
[873]370  return
371
372
373end subroutine optcv
Note: See TracBrowser for help on using the repository browser.