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

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