source: trunk/LMDZ.TITAN/libf/phytitan/optcv.F90 @ 1722

Last change on this file since 1722 was 1722, checked in by jvatant, 7 years ago

Adapt various modifs of LMDZ.GENERIC to LMDZ.TITAN from r1690-1694-1699-1709-1715
--JVO

  • Property svn:executable set to *
File size: 11.9 KB
Line 
1SUBROUTINE OPTCV(DTAUV,TAUV,TAUCUMV,PLEV,  &
2     QXVAER,QSVAER,GVAER,WBARV,COSBV,       &
3     TAURAY,TAUAERO,TMID,PMID,TAUGSURF,GWEIGHT)
4
5  use radinc_h
6  use radcommon_h, only: gasv, tlimit, Cmk, tgasref, pfgasref,wnov,scalep,indv,glat_ig
7  use gases_h
8  use comcstfi_mod, only: g, r
9  use callkeys_mod, only: continuum,graybody,callgasvis
10
11  implicit none
12
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 CALCULATES FOR EACH LAYER, FOR EACH SPECTRAL 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  !-------------------------------------------------------------------
38
39
40  real*8 DTAUV(L_NLAYRAD,L_NSPECTV,L_NGAUSS)
41  real*8 DTAUKV(L_LEVELS,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)
48
49  ! for aerosols
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)
55  real*8  TAEROS(L_LEVELS,L_NSPECTV,NAERKIND)
56
57  ! Titan customisation
58  ! J. Vatant d'Ollone (2016)
59  real*8 GWEIGHT(L_NGAUSS)
60  real*8 DHAZE_T(L_LEVELS,L_NSPECTI)
61  real*8 DHAZES_T(L_LEVELS,L_NSPECTI)
62  real*8 SSA_T(L_LEVELS,L_NSPECTI)
63  real*8 ASF_T(L_LEVELS,L_NSPECTI)
64  real*8 INT_DTAU(L_NLAYRAD,L_NSPECTI)
65  real*8 K_HAZE(L_NLAYRAD,L_NSPECTI)
66 
67  CHARACTER*2  str2
68  ! ==========================
69
70  integer L, NW, NG, K, LK, IAER
71  integer MT(L_LEVELS), MP(L_LEVELS), NP(L_LEVELS)
72  real*8  ANS, TAUGAS
73  real*8  TAURAY(L_NSPECTV)
74  real*8  TRAY(L_LEVELS,L_NSPECTV)
75  real*8  DPR(L_LEVELS), U(L_LEVELS)
76  real*8  LCOEF(4), LKCOEF(L_LEVELS,4)
77
78  real*8 taugsurf(L_NSPECTV,L_NGAUSS-1)
79  real*8 DCONT,DAERO
80  real*8 DRAYAER
81  double precision wn_cont, p_cont, p_air, T_cont, dtemp, dtempc
82  double precision p_cross
83
84  real*8  KCOEF(4)
85
86  ! temporary variables for multiple aerosol calculation
87  real*8 atemp(L_NLAYRAD,L_NSPECTV)
88  real*8 btemp(L_NLAYRAD,L_NSPECTV)
89  real*8 ctemp(L_NLAYRAD,L_NSPECTV)
90
91  ! variables for k in units m^-1
92  real*8 dz(L_LEVELS)
93
94  integer igas, jgas, ilay
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
114  do K=2,L_LEVELS
115     DPR(k) = PLEV(K)-PLEV(K-1)
116
117     ! if we have continuum opacities, we need dz
118
119      dz(k) = dpr(k)*R*TMID(K)/(glat_ig*PMID(K))
120      U(k)  = Cmk*DPR(k)     ! only Cmk line in optcv.F     
121
122     call tpindex(PMID(K),TMID(K),pfgasref,tgasref,LCOEF,MT(K),MP(K))
123
124     do LK=1,4
125        LKCOEF(K,LK) = LCOEF(LK)
126     end do
127  end do                    ! levels
128
129  ! Spectral dependance of aerosol absorption
130  do iaer=1,naerkind
131     do NW=1,L_NSPECTV
132        do K=2,L_LEVELS
133           TAEROS(K,NW,IAER) = TAUAERO(K,IAER) * QXVAER(K,NW,IAER)
134        end do                    ! levels
135     end do
136  end do
137 
138  ! Rayleigh scattering
139  do NW=1,L_NSPECTV
140     do K=2,L_LEVELS
141        TRAY(K,NW)   = TAURAY(NW) * DPR(K)
142     end do                    ! levels
143  end do
144 
145  !     we ignore K=1...
146  do K=2,L_LEVELS
147 
148     ilay = k / 2 ! int. arithmetic => gives the gcm layer index
149
150     do NW=1,L_NSPECTV
151
152        !================= Titan customisation ========================================
153        call disr_haze(dz(k),plev(k),wnov(nw),dhaze_T(k,nw),SSA_T(k,nw),ASF_T(k,nw))
154        ! =============================================================================
155
156        DRAYAER = TRAY(K,NW)
157        !     DRAYAER is Tau RAYleigh scattering, plus AERosol opacity
158        do iaer=1,naerkind
159           DRAYAER = DRAYAER + TAEROS(K,NW,IAER)
160        end do
161
162        DRAYAER = DRAYAER + DHAZE_T(K,NW) ! Titan's aerosol
163
164        DCONT = 0.0 ! continuum absorption
165
166        if(continuum.and.(.not.graybody).and.callgasvis)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              p_cont  = dble(PMID(k)*scalep*gfrac(igas,ilay))
173
174              dtemp=0.0
175              if(igas.eq.igas_N2)then
176
177                 interm = indv(nw,igas,igas)
178!                 call interpolateN2N2(wn_cont,T_cont,p_cont,dtemp,.false.,interm)
179                 indv(nw,igas,igas) = interm
180                 ! only goes to 500 cm^-1, so unless we're around a cold brown dwarf, this is irrelevant in the visible
181
182              elseif(igas.eq.igas_H2)then
183
184                 ! first do self-induced absorption
185                 interm = indv(nw,igas,igas)
186                 call interpolateH2H2(wn_cont,T_cont,p_cont,dtemp,.false.,interm)
187                 indv(nw,igas,igas) = interm
188
189                 ! then cross-interactions with other gases
190                 do jgas=1,ngasmx
191                    p_cross = dble(PMID(k)*scalep*gfrac(jgas,ilay))
192                    dtempc  = 0.0
193                    if(jgas.eq.igas_N2)then
194                       interm = indv(nw,igas,jgas)
195                       call interpolateN2H2(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm)
196                       indv(nw,igas,jgas) = interm
197                       ! should be irrelevant in the visible
198                    endif
199                    dtemp = dtemp + dtempc
200                 enddo
201
202               elseif(igas.eq.igas_CH4)then
203
204                 ! first do self-induced absorption
205                 interm = indv(nw,igas,igas)
206                 call interpolateCH4CH4(wn_cont,T_cont,p_cont,dtemp,.false.,interm)
207                 indv(nw,igas,igas) = interm
208
209                 ! then cross-interactions with other gases
210                 do jgas=1,ngasmx
211                    p_cross = dble(PMID(k)*scalep*gfrac(jgas,ilay))
212                    dtempc  = 0.0
213                    if(jgas.eq.igas_N2)then
214                       interm = indv(nw,igas,jgas)
215                       call interpolateN2CH4(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm)
216                       indv(nw,igas,jgas) = interm
217                    endif
218                    dtemp = dtemp + dtempc
219                 enddo
220
221              endif
222
223              DCONT = DCONT + dtemp
224
225           enddo
226
227           DCONT = DCONT*dz(k)
228
229        endif
230
231        do ng=1,L_NGAUSS-1
232
233           ! Now compute TAUGAS
234
235           KCOEF(1) = GASV(MT(K),MP(K),1,NW,NG)
236           KCOEF(2) = GASV(MT(K),MP(K)+1,1,NW,NG)
237           KCOEF(3) = GASV(MT(K)+1,MP(K)+1,1,NW,NG)
238           KCOEF(4) = GASV(MT(K)+1,MP(K),1,NW,NG)
239
240           ! Interpolate the gaseous k-coefficients to the requested T,P values
241
242           ANS = LKCOEF(K,1)*KCOEF(1) + LKCOEF(K,2)*KCOEF(2) +            &
243                LKCOEF(K,3)*KCOEF(3) + LKCOEF(K,4)*KCOEF(4)
244
245
246           TAUGAS  = U(k)*ANS
247
248           TAUGSURF(NW,NG) = TAUGSURF(NW,NG) + TAUGAS + DCONT
249           DTAUKV(K,nw,ng) = TAUGAS &
250                             + DRAYAER & ! DRAYAER includes all scattering contributions
251                             + DCONT ! For parameterized continuum aborption
252
253        end do
254
255        ! Now fill in the "clear" part of the spectrum (NG = L_NGAUSS),
256        ! which holds continuum opacity only
257
258        NG              = L_NGAUSS
259        DTAUKV(K,nw,ng) = DRAYAER + DCONT ! Scattering + parameterized continuum absorption, including Titan's haze
260
261     end do
262  end do
263
264
265  !=======================================================================
266  !     Now the full treatment for the layers, where besides the opacity
267  !     we need to calculate the scattering albedo and asymmetry factors
268
269  do iaer=1,naerkind
270    DO NW=1,L_NSPECTV
271      DO K=2,L_LEVELS
272           TAUAEROLK(K,NW,IAER) = TAUAERO(K,IAER) * QSVAER(K,NW,IAER) ! effect of scattering albedo
273      ENDDO
274    ENDDO
275  end do
276 
277  ! Haze scattering
278  DO NW=1,L_NSPECTV
279    DO K=2,L_LEVELS
280      DHAZES_T(K,NW) = DHAZE_T(K,NW) * SSA_T(K,NW) ! effect of scattering albedo on haze
281    ENDDO
282  ENDDO
283
284
285  DO NW=1,L_NSPECTV
286     DO L=1,L_NLAYRAD-1
287        K              = 2*L+1
288        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)) &
289                    + ASF_T(K,NW)*DHAZES_T(K,NW) + ASF_T(K+1,NW)*DHAZES_T(K+1,NW)
290        btemp(L,NW) = SUM(TAUAEROLK(K,NW,1:naerkind)) + SUM(TAUAEROLK(K+1,NW,1:naerkind)) &
291                    + DHAZES_T(K,NW) + DHAZES_T(K+1,NW)
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                 + ASF_T(K,NW)*DHAZES_T(K,NW)
302     btemp(L,NW) = SUM(TAUAEROLK(K,NW,1:naerkind)) &
303                 + DHAZES_T(K,NW)
304     ctemp(L,NW) = btemp(L,NW) + 0.9999*TRAY(K,NW) ! JVO 2017 : does this 0.999 is really meaningful ?
305     btemp(L,NW) = btemp(L,NW) + TRAY(K,NW)
306     COSBV(L,NW,1:L_NGAUSS) = atemp(L,NW)/btemp(L,NW)
307     
308     
309  END DO                    ! NW spectral loop
310
311  DO NG=1,L_NGAUSS
312    DO NW=1,L_NSPECTV
313     DO L=1,L_NLAYRAD-1
314
315        K              = 2*L+1
316        DTAUV(L,nw,ng) = DTAUKV(K,NW,NG) + DTAUKV(K+1,NW,NG)
317        WBARV(L,nw,ng) = ctemp(L,NW) / DTAUV(L,nw,ng)
318
319      END DO ! L vertical loop
320
321        ! Last level
322
323        L              = L_NLAYRAD
324        K              = 2*L+1
325        DTAUV(L,nw,ng) = DTAUKV(K,NW,NG)
326
327        WBARV(L,NW,NG) = ctemp(L,NW) / DTAUV(L,NW,NG)
328
329     END DO                 ! NW spectral loop
330  END DO                    ! NG Gauss loop
331
332  ! Total extinction optical depths
333
334  DO NG=1,L_NGAUSS       ! full gauss loop
335     DO NW=1,L_NSPECTV       
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
345     END DO           
346  END DO                 ! end full gauss loop
347
348
349!  Titan's outputs (JVO, 2016)===============================================
350!      do l=1,L_NLAYRAD
351!         do nw=1,L_NSPECTV
352!          INT_DTAU(L,NW) = 0.0d+0
353!            DO NG=1,L_NGAUSS
354!               INT_DTAU(L,NW)= INT_DTAU(L,NW) + dtauv(L,nw,ng)*gweight(NG)
355!            enddo
356!         enddo
357!      enddo
358
359!       do nw=1,L_NSPECTV
360!          write(str2,'(i2.2)') nw
361!         call writediagfi(1,'kgv'//str2,'Gaz extinction coefficient VI band '//str2,'m-1',1,int_dtau(L_NLAYRAD:1:-1,nw)/dz_lay(L_NLAYRAD:1:-1))
362!          call writediagfi(1,'khv'//str2,'Haze extinction coefficient VI band '//str2,'m-1',1,k_haze(L_NLAYRAD:1:-1,nw)/dz_lay(L_NLAYRAD:1:-1))       
363!       enddo 
364
365! ============================================================================== 
366
367
368  return
369
370
371end subroutine optcv
Note: See TracBrowser for help on using the repository browser.