source: trunk/LMDZ.VENUS/libf/phyvenus/optcv.F90 @ 3556

Last change on this file since 3556 was 2560, checked in by slebonnois, 3 years ago

SL: Implementation of SW computation based on generic model. Switch between this new SW module or old module that reads R. Haus tables implemented with a key (solarchoice)

File size: 11.2 KB
Line 
1MODULE optcv_mod
2
3IMPLICIT NONE
4
5CONTAINS
6
7SUBROUTINE OPTCV(DTAUV,TAUV,TAUCUMV,PLEV,  &
8     QXVAER,QSVAER,GVAER,WBARV,COSBV,       &
9     TAURAY,TAUAERO,TMID,PMID,TAUGSURF,QVAR,MUVAR)
10
11  use radinc_h, only: L_NLAYRAD, L_NLEVRAD, L_LEVELS, L_NSPECTV, L_NGAUSS, L_REFVAR, NAERKIND
12  use radcommon_h, only: gasv, tlimit, tgasref, pfgasref,wnov,scalep,indv
13  use gases_h, only: gfrac, ngasmx, igas_H2, igas_H2O, igas_He, igas_N2
14
15  implicit none
16#include "YOMCST.h"
17
18  !==================================================================
19  !     
20  !     Purpose
21  !     -------
22  !     Calculates shortwave optical constants at each level.
23  !     
24  !     Authors
25  !     -------
26  !     Adapted from the NASA Ames code by R. Wordsworth (2009)
27  !     
28  !==================================================================
29  !     
30  !     THIS SUBROUTINE SETS THE OPTICAL CONSTANTS IN THE VISUAL 
31  !     IT CALCULATES FOR EACH LAYER, FOR EACH SPECTRAL INTERVAL IN THE VISUAL
32  !     LAYER: WBAR, DTAU, COSBAR
33  !     LEVEL: TAU
34  !     
35  !     TAUV(L,NW,NG) is the cumulative optical depth at the top of radiation code
36  !     layer L. NW is spectral wavelength interval, ng the Gauss point index.
37  !     
38  !     TLEV(L) - Temperature at the layer boundary
39  !     PLEV(L) - Pressure at the layer boundary (i.e. level)
40  !     GASV(NT,NPS,NW,NG) - Visible k-coefficients
41  !     
42  !-------------------------------------------------------------------
43
44
45  real*8 DTAUV(L_NLAYRAD,L_NSPECTV,L_NGAUSS)
46  real*8 DTAUKV(L_LEVELS,L_NSPECTV,L_NGAUSS)
47  real*8 TAUV(L_NLEVRAD,L_NSPECTV,L_NGAUSS)
48  real*8 TAUCUMV(L_LEVELS,L_NSPECTV,L_NGAUSS)
49  real*8 PLEV(L_LEVELS)
50  real*8 TMID(L_LEVELS), PMID(L_LEVELS)
51  real*8 COSBV(L_NLAYRAD,L_NSPECTV,L_NGAUSS)
52  real*8 WBARV(L_NLAYRAD,L_NSPECTV,L_NGAUSS)
53
54  ! for aerosols
55  real*8  QXVAER(L_LEVELS,L_NSPECTV,NAERKIND)
56  real*8  QSVAER(L_LEVELS,L_NSPECTV,NAERKIND)
57  real*8  GVAER(L_LEVELS,L_NSPECTV,NAERKIND)
58  real*8  TAUAERO(L_LEVELS,NAERKIND)
59  real*8  TAUAEROLK(L_LEVELS,L_NSPECTV,NAERKIND)
60  real*8  TAEROS(L_LEVELS,L_NSPECTV,NAERKIND)
61
62  integer L, NW, NG, K, LK, IAER
63  integer MT(L_LEVELS), MP(L_LEVELS), NP(L_LEVELS)
64  real*8  ANS, TAUGAS
65  real*8  TAURAY(L_NSPECTV)
66  real*8  TRAY(L_LEVELS,L_NSPECTV)
67  real*8  DPR(L_LEVELS), U(L_LEVELS)
68  real*8  LCOEF(4), LKCOEF(L_LEVELS,4)
69
70  real*8 taugsurf(L_NSPECTV,L_NGAUSS-1)
71  real*8 DCONT,DAERO
72  real*8 DRAYAER
73  double precision wn_cont, p_cont, p_air, T_cont, dtemp, dtempc
74  double precision p_cross
75
76  ! variable species mixing ratio variables
77  real*8  QVAR(L_LEVELS), WRATIO(L_LEVELS), MUVAR(L_LEVELS)
78  real*8  KCOEF(4)
79  integer NVAR(L_LEVELS)
80 
81  ! temporary variables to reduce memory access time to gasv
82  real*8 tmpk(2,2)
83  real*8 tmpkvar(2,2,2)
84
85  ! temporary variables for multiple aerosol calculation
86  real*8 atemp(L_NLAYRAD,L_NSPECTV)
87  real*8 btemp(L_NLAYRAD,L_NSPECTV)
88  real*8 ctemp(L_NLAYRAD,L_NSPECTV)
89
90  ! variables for k in units m^-1
91  real*8 dz(L_LEVELS)
92  real*8 Cmk
93
94  integer igas, jgas
95
96  integer interm
97
98  !! AS: to save time in computing continuum (see bilinearbig)
99  IF (.not.ALLOCATED(indv)) THEN
100      ALLOCATE(indv(L_NSPECTV,ngasmx,ngasmx))
101      indv = -9999 ! this initial value means "to be calculated"
102  ENDIF
103
104  !=======================================================================
105  !     Determine the total gas opacity throughout the column, for each
106  !     spectral interval, NW, and each Gauss point, NG.
107  !     Calculate the continuum opacities, i.e., those that do not depend on
108  !     NG, the Gauss index.
109
110  taugsurf(:,:) = 0.0
111  dpr(:)        = 0.0
112  lkcoef(:,:)   = 0.0
113  Cmk= 0.01 * 1.0 / (RG * RMD * 1.672621e-27) ! q_main=1.0 assumed.
114
115  do K=2,L_LEVELS
116     DPR(k) = PLEV(K)-PLEV(K-1)
117
118     ! if we have continuum opacities, we need dz
119        dz(k) = dpr(k)*R*TMID(K)/(RG*PMID(K))
120        U(k)  = Cmk*DPR(k)
121
122     call tpindex(PMID(K),TMID(K),QVAR(K),pfgasref,tgasref,&
123          LCOEF,MT(K),MP(K),NVAR(K),WRATIO(K))
124
125     do LK=1,4
126        LKCOEF(K,LK) = LCOEF(LK)
127     end do
128  end do                    ! levels
129
130  ! Spectral dependance of aerosol absorption
131            !JL18 It seems to be good to have aerosols in the first "radiative layer" of the gcm in the IR
132            !   but visible does not handle very well diffusion in first layer.
133            !   The tauaero and tauray are thus set to 0 (a small value for rayleigh because the code crashes otherwise)
134            !   in the 4 first semilayers in optcv, but not optci.
135            !   This solves random variations of the sw heating at the model top.
136  do iaer=1,naerkind
137     do NW=1,L_NSPECTV
138        TAEROS(1:4,NW,IAER)=0.d0
139        do K=5,L_LEVELS
140           TAEROS(K,NW,IAER) = TAUAERO(K,IAER) * QXVAER(K,NW,IAER)
141        end do                    ! levels
142     end do
143  end do
144 
145  ! Rayleigh scattering
146  do NW=1,L_NSPECTV
147     TRAY(1:4,NW)   = 1d-30
148     do K=5,L_LEVELS
149        TRAY(K,NW)   = TAURAY(NW) * DPR(K)
150     end do                    ! levels
151  end do
152 
153  !     we ignore K=1...
154  do K=2,L_LEVELS
155
156     do NW=1,L_NSPECTV
157
158        DRAYAER = TRAY(K,NW)
159        !     DRAYAER is Tau RAYleigh scattering, plus AERosol opacity
160        do iaer=1,naerkind
161           DRAYAER = DRAYAER + TAEROS(K,NW,IAER)
162        end do
163
164        DCONT = 0.0 ! continuum absorption
165
166!        if(continuum)then
167           ! include continua if necessary
168           wn_cont = dble(wnov(nw))
169           T_cont  = dble(TMID(k))
170           do igas=1,ngasmx
171
172              if(gfrac(igas).eq.-1)then ! variable
173                 p_cont  = dble(PMID(k)*scalep*QVAR(k)) ! qvar = mol/mol
174              else
175                 p_cont  = dble(PMID(k)*scalep*gfrac(igas)*(1.-QVAR(k)))
176              endif
177
178              dtemp=0.0
179
180! For Venus: only H2O, using CKD
181              if(igas.eq.igas_H2O.and.T_cont.gt.200.0)then
182
183                 p_air = dble(PMID(k)*scalep) - p_cont ! note assumes background is air!
184!                 if(H2Ocont_simple)then
185!                    call interpolateH2Ocont_PPC(wn_cont,T_cont,p_cont,p_air,dtemp,.false.)
186!                 else
187                    interm = indv(nw,igas,igas)
188                    call interpolateH2Ocont_CKD(wn_cont,T_cont,p_cont,p_air,dtemp,.false.,interm)
189                    indv(nw,igas,igas) = interm
190!                 endif
191
192              endif
193
194              DCONT = DCONT + dtemp
195
196           enddo
197
198           DCONT = DCONT*dz(k)
199
200!        endif
201
202        do ng=1,L_NGAUSS-1
203
204           ! Now compute TAUGAS
205
206           ! Interpolate between water mixing ratios
207           ! WRATIO = 0.0 if the requested water amount is equal to, or outside the
208           ! the water data range
209
210           if(L_REFVAR.eq.1)then ! added by RW for special no variable case
211           
212              ! JVO 2017 : added tmpk because the repeated calls to gasi/v increased dramatically
213              ! the execution time of optci/v -> ~ factor 2 on the whole radiative
214              ! transfer on the tested simulations !
215
216              tmpk = GASV(MT(K):MT(K)+1,MP(K):MP(K)+1,1,NW,NG)
217             
218              KCOEF(1) = tmpk(1,1) ! KCOEF(1) = GASV(MT(K),MP(K),1,NW,NG)
219              KCOEF(2) = tmpk(1,2) ! KCOEF(2) = GASV(MT(K),MP(K)+1,1,NW,NG)
220              KCOEF(3) = tmpk(2,2) ! KCOEF(3) = GASV(MT(K)+1,MP(K)+1,1,NW,NG)
221              KCOEF(4) = tmpk(2,1) ! KCOEF(4) = GASV(MT(K)+1,MP(K),1,NW,NG)
222
223           else
224
225              tmpkvar = GASV(MT(K):MT(K)+1,MP(K):MP(K)+1,NVAR(K):NVAR(K)+1,NW,NG)
226
227              KCOEF(1) = tmpkvar(1,1,1) + WRATIO(K) *  &
228                        ( tmpkvar(1,1,2)-tmpkvar(1,1,1) )
229
230              KCOEF(2) = tmpkvar(1,2,1) + WRATIO(K) *  &
231                        ( tmpkvar(1,2,2)-tmpkvar(1,2,1) )
232
233              KCOEF(3) = tmpkvar(2,2,1) + WRATIO(K) *  &
234                        ( tmpkvar(2,2,2)-tmpkvar(2,2,1) )
235             
236              KCOEF(4) = tmpkvar(2,1,1) + WRATIO(K) *  &
237                        ( tmpkvar(2,1,2)-tmpkvar(2,1,1) )
238
239
240           endif
241
242           ! Interpolate the gaseous k-coefficients to the requested T,P values
243
244           ANS = LKCOEF(K,1)*KCOEF(1) + LKCOEF(K,2)*KCOEF(2) +            &
245                LKCOEF(K,3)*KCOEF(3) + LKCOEF(K,4)*KCOEF(4)
246
247           TAUGAS  = U(k)*ANS
248
249           TAUGSURF(NW,NG) = TAUGSURF(NW,NG) + TAUGAS + DCONT
250           DTAUKV(K,nw,ng) = TAUGAS &
251                             + DRAYAER & ! DRAYAER includes all scattering contributions
252                             + DCONT ! For parameterized continuum aborption
253
254        end do
255
256        ! Now fill in the "clear" part of the spectrum (NG = L_NGAUSS),
257        ! which holds continuum opacity only
258
259        NG              = L_NGAUSS
260        DTAUKV(K,nw,ng) = DRAYAER + DCONT ! Scattering + parameterized continuum absorption
261
262     end do
263  end do
264
265
266  !=======================================================================
267  !     Now the full treatment for the layers, where besides the opacity
268  !     we need to calculate the scattering albedo and asymmetry factors
269
270            !JL18 It seems to be good to have aerosols in the first "radiative layer" of the gcm in the IR
271            !   but not in the visible
272            !   The tauaero is thus set to 0 in the 4 first semilayers in optcv, but not optci.
273            !   This solves random variations of the sw heating at the model top.
274  do iaer=1,naerkind
275    DO NW=1,L_NSPECTV
276      TAUAEROLK(1:4,NW,IAER)=0.d0
277      DO K=5,L_LEVELS
278           TAUAEROLK(K,NW,IAER) = TAUAERO(K,IAER) * QSVAER(K,NW,IAER) ! effect of scattering albedo
279      ENDDO
280    ENDDO
281  end do
282
283  DO NW=1,L_NSPECTV
284     DO L=1,L_NLAYRAD-1
285        K              = 2*L+1
286        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))
287        btemp(L,NW) = SUM(TAUAEROLK(K,NW,1:naerkind)) + SUM(TAUAEROLK(K+1,NW,1:naerkind))
288        ctemp(L,NW) = btemp(L,NW) + 0.9999*(TRAY(K,NW) + TRAY(K+1,NW))  ! JVO 2017 : does this 0.999 is really meaningful ?
289        btemp(L,NW) = btemp(L,NW) + TRAY(K,NW) + TRAY(K+1,NW)
290        COSBV(L,NW,1:L_NGAUSS) = atemp(L,NW)/btemp(L,NW)
291     END DO ! L vertical loop
292     
293     ! Last level
294     L           = L_NLAYRAD
295     K           = 2*L+1
296     atemp(L,NW) = SUM(GVAER(K,NW,1:naerkind) * TAUAEROLK(K,NW,1:naerkind))
297     btemp(L,NW) = SUM(TAUAEROLK(K,NW,1:naerkind))
298     ctemp(L,NW) = btemp(L,NW) + 0.9999*TRAY(K,NW) ! JVO 2017 : does this 0.999 is really meaningful ?
299     btemp(L,NW) = btemp(L,NW) + TRAY(K,NW)
300     COSBV(L,NW,1:L_NGAUSS) = atemp(L,NW)/btemp(L,NW)
301     
302     
303  END DO                    ! NW spectral loop
304
305  DO NG=1,L_NGAUSS
306    DO NW=1,L_NSPECTV
307     DO L=1,L_NLAYRAD-1
308
309        K              = 2*L+1
310        DTAUV(L,nw,ng) = DTAUKV(K,NW,NG) + DTAUKV(K+1,NW,NG)
311        WBARV(L,nw,ng) = ctemp(L,NW) / DTAUV(L,nw,ng)
312
313      END DO ! L vertical loop
314
315        ! Last level
316
317        L              = L_NLAYRAD
318        K              = 2*L+1
319        DTAUV(L,nw,ng) = DTAUKV(K,NW,NG)
320
321        WBARV(L,NW,NG) = ctemp(L,NW) / DTAUV(L,NW,NG)
322
323     END DO                 ! NW spectral loop
324  END DO                    ! NG Gauss loop
325
326  ! Total extinction optical depths
327
328  DO NG=1,L_NGAUSS       ! full gauss loop
329     DO NW=1,L_NSPECTV       
330        TAUCUMV(1,NW,NG)=0.0D0
331        DO K=2,L_LEVELS
332           TAUCUMV(K,NW,NG)=TAUCUMV(K-1,NW,NG)+DTAUKV(K,NW,NG)
333        END DO
334
335        DO L=1,L_NLAYRAD
336           TAUV(L,NW,NG)=TAUCUMV(2*L,NW,NG)
337        END DO
338        TAUV(L,NW,NG)=TAUCUMV(2*L_NLAYRAD+1,NW,NG)
339     END DO           
340  END DO                 ! end full gauss loop
341
342
343
344
345end subroutine optcv
346
347END MODULE optcv_mod
Note: See TracBrowser for help on using the repository browser.