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

Last change on this file since 1648 was 1648, checked in by jvatant, 8 years ago

Modifications to custom radiative transfer to Titan
+ Enables an altitude dependant gfrac for CIA computations

-> many radical changes in su_gases and co ..
-> read vertical CH4 profile with call_profilgases
-> Now you need a 'profile.def' that I will add in the deftank

+ Added interpolate CIA routines for CH4
+ Added temporary mean aerosol profile opacity routine (disr_haze)

  • Property svn:executable set to *
File size: 12.1 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 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  !-------------------------------------------------------------------
38
39
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)
48
49  ! for aerosols
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)
56
57  ! Titan customisation
58  ! J. Vatant d'Ollone (2016)
59  real*8 GWEIGHT(L_NGAUSS)
60  real*8 DHAZE_T(L_LEVELS+1,L_NSPECTI)
61  real*8 DHAZES_T(L_LEVELS+1,L_NSPECTI)
62  real*8 SSA_T(L_LEVELS+1,L_NSPECTI)
63  real*8 ASF_T(L_LEVELS+1,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  TRAYAER
76  real*8  DPR(L_LEVELS), U(L_LEVELS)
77  real*8  LCOEF(4), LKCOEF(L_LEVELS,4)
78
79  real*8 taugsurf(L_NSPECTV,L_NGAUSS-1)
80  real*8 DCONT,DAERO
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
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  do NW=1,L_NSPECTV
138     do K=2,L_LEVELS
139        TRAY(K,NW)   = TAURAY(NW) * DPR(K)
140     end do                    ! levels
141  end do
142
143  !     we ignore K=1...
144  do K=2,L_LEVELS
145 
146     ilay = k / 2 ! int. arithmetic => gives the gcm layer index
147
148     do NW=1,L_NSPECTV
149
150        TRAYAER = TRAY(K,NW)
151        !     TRAYAER is Tau RAYleigh scattering, plus AERosol opacity
152        do iaer=1,naerkind
153           TRAYAER = TRAYAER + TAEROS(K,NW,IAER)
154        end do
155
156        DCONT = 0.0 ! continuum absorption
157
158        if(continuum.and.(.not.graybody).and.callgasvis)then
159           ! include continua if necessary
160           wn_cont = dble(wnov(nw))
161           T_cont  = dble(TMID(k))
162           do igas=1,ngasmx
163
164              p_cont  = dble(PMID(k)*scalep*gfrac(igas,ilay))
165
166              dtemp=0.0
167              if(igas.eq.igas_N2)then
168
169                 interm = indv(nw,igas,igas)
170!                 call interpolateN2N2(wn_cont,T_cont,p_cont,dtemp,.false.,interm)
171                 indv(nw,igas,igas) = interm
172                 ! only goes to 500 cm^-1, so unless we're around a cold brown dwarf, this is irrelevant in the visible
173
174              elseif(igas.eq.igas_H2)then
175
176                 ! first do self-induced absorption
177                 interm = indv(nw,igas,igas)
178                 call interpolateH2H2(wn_cont,T_cont,p_cont,dtemp,.false.,interm)
179                 indv(nw,igas,igas) = interm
180
181                 ! then cross-interactions with other gases
182                 do jgas=1,ngasmx
183                    p_cross = dble(PMID(k)*scalep*gfrac(jgas,ilay))
184                    dtempc  = 0.0
185                    if(jgas.eq.igas_N2)then
186                       interm = indv(nw,igas,jgas)
187                       call interpolateN2H2(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm)
188                       indv(nw,igas,jgas) = interm
189                       ! should be irrelevant in the visible
190                    endif
191                    dtemp = dtemp + dtempc
192                 enddo
193
194               elseif(igas.eq.igas_CH4)then
195
196                 ! first do self-induced absorption
197                 interm = indv(nw,igas,igas)
198                 call interpolateCH4CH4(wn_cont,T_cont,p_cont,dtemp,.false.,interm)
199                 indv(nw,igas,igas) = interm
200
201                 ! then cross-interactions with other gases
202                 do jgas=1,ngasmx
203                    p_cross = dble(PMID(k)*scalep*gfrac(jgas,ilay))
204                    dtempc  = 0.0
205                    if(jgas.eq.igas_N2)then
206                       interm = indv(nw,igas,jgas)
207                       call interpolateN2CH4(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm)
208                       indv(nw,igas,jgas) = interm
209                    endif
210                    dtemp = dtemp + dtempc
211                 enddo
212
213              endif
214
215              DCONT = DCONT + dtemp
216
217           enddo
218
219           DCONT = DCONT*dz(k)
220
221        endif
222
223!================= Titan customisation ========================================
224        call disr_haze(dz(k),plev(k),wnov(nw),dhaze_T(k,nw),SSA_T(k,nw),ASF_T(k,nw))
225! =============================================================================
226
227        do ng=1,L_NGAUSS-1
228
229           ! Now compute TAUGAS
230
231           KCOEF(1) = GASV(MT(K),MP(K),1,NW,NG)
232           KCOEF(2) = GASV(MT(K),MP(K)+1,1,NW,NG)
233           KCOEF(3) = GASV(MT(K)+1,MP(K)+1,1,NW,NG)
234           KCOEF(4) = GASV(MT(K)+1,MP(K),1,NW,NG)
235
236           ! Interpolate the gaseous k-coefficients to the requested T,P values
237
238           ANS = LKCOEF(K,1)*KCOEF(1) + LKCOEF(K,2)*KCOEF(2) +            &
239                LKCOEF(K,3)*KCOEF(3) + LKCOEF(K,4)*KCOEF(4)
240
241           TAUGAS  = U(k)*ANS
242
243           TAUGSURF(NW,NG) = TAUGSURF(NW,NG) + TAUGAS + DCONT
244           DTAUKV(K,nw,ng) = TAUGAS &
245                             + TRAYAER & ! TRAYAER includes all scattering contributions
246                             + DCONT    & ! For parameterized continuum aborption
247                             + DHAZE_T(K,NW)  ! For Titan haze
248
249        end do
250
251        ! Now fill in the "clear" part of the spectrum (NG = L_NGAUSS),
252        ! which holds continuum opacity only
253
254        NG              = L_NGAUSS
255        DTAUKV(K,nw,ng) = TRAY(K,NW) + DCONT  & ! For parameterized continuum absorption
256                              + DHAZE_T(K,NW)  ! For Titan haze
257
258        do iaer=1,naerkind
259           DTAUKV(K,nw,ng) = DTAUKV(K,nw,ng) +  TAEROS(K,NW,IAER)
260        end do ! a bug was here!
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 
271  do iaer=1,naerkind
272    DO NW=1,L_NSPECTV
273      DO K=2,L_LEVELS   ! AS: shouldn't this be L_LEVELS+1 ? (see optci)
274           TAUAEROLK(K,NW,IAER) = TAUAERO(K,IAER) * QSVAER(K,NW,IAER)
275      ENDDO
276    ENDDO
277  end do
278 
279  ! Haze scattering
280  DO NW=1,L_NSPECTV
281    DO K=2,L_LEVELS+1
282      DHAZES_T(K,NW) = DHAZE_T(K,NW) * SSA_T(K,NW)
283    ENDDO
284  ENDDO
285
286
287  DO NW=1,L_NSPECTV
288     DO L=1,L_NLAYRAD-1
289        K              = 2*L+1
290       
291        atemp(L,NW) = SUM(GVAER(K,NW,1:naerkind) * TAUAEROLK(K,NW,1:naerkind)) &
292                    + SUM(GVAER(K+1,NW,1:naerkind) * TAUAEROLK(K+1,NW,1:naerkind)) &
293                    + ASF_T(K,NW)*DHAZES_T(K,NW)
294                   
295        btemp(L,NW) = SUM(TAUAEROLK(K,NW,1:naerkind)) + SUM(TAUAEROLK(K+1,NW,1:naerkind)) &
296                    + DHAZES_T(K,NW)
297       
298        ctemp(L,NW) = btemp(L,NW) + 0.9999*(TRAY(K,NW) + TRAY(K+1,NW))
299        btemp(L,NW) = btemp(L,NW) + TRAY(K,NW) + TRAY(K+1,NW)
300        COSBV(L,NW,1:L_NGAUSS) = atemp(L,NW)/btemp(L,NW)
301     END DO ! L vertical loop
302     
303     !last level
304     L              = L_NLAYRAD
305     K              = 2*L+1
306     
307     atemp(L,NW) = SUM(GVAER(K,NW,1:naerkind) * TAUAEROLK(K,NW,1:naerkind)) &
308                 + ASF_T(K,NW)*DHAZES_T(K,NW)
309     btemp(L,NW) = SUM(TAUAEROLK(K,NW,1:naerkind)) &
310                 + DHAZES_T(K,NW)
311     ctemp(L,NW) = btemp(L,NW) + 0.9999*TRAY(K,NW)
312     btemp(L,NW) = btemp(L,NW) + TRAY(K,NW)
313     COSBV(L,NW,1:L_NGAUSS) = atemp(L,NW)/btemp(L,NW)
314     
315     
316  END DO                    ! NW spectral loop
317
318! ===========================================================================================
319
320  DO NG=1,L_NGAUSS
321    DO NW=1,L_NSPECTV
322     DO L=1,L_NLAYRAD-1
323
324        K              = 2*L+1
325        DTAUV(L,nw,ng) = DTAUKV(K,NW,NG) + DTAUKV(K+1,NW,NG)
326        WBARV(L,nw,ng) = ctemp(L,NW) / DTAUV(L,nw,ng)
327
328      END DO ! L vertical loop
329
330        !     No vertical averaging on bottom layer
331
332        L              = L_NLAYRAD
333        K              = 2*L+1
334        DTAUV(L,nw,ng) = DTAUKV(K,NW,NG)
335
336        WBARV(L,NW,NG) = ctemp(L,NW) / DTAUV(L,NW,NG)
337
338     END DO                 ! NW spectral loop
339  END DO                    ! NG Gauss loop
340
341  ! Total extinction optical depths
342
343  DO NG=1,L_NGAUSS       ! full gauss loop
344     DO NW=1,L_NSPECTV       
345        TAUV(1,NW,NG)=0.0D0
346        DO L=1,L_NLAYRAD
347           TAUV(L+1,NW,NG)=TAUV(L,NW,NG)+DTAUV(L,NW,NG)
348        END DO
349
350        TAUCUMV(1,NW,NG)=0.0D0
351        DO K=2,L_LEVELS
352           TAUCUMV(K,NW,NG)=TAUCUMV(K-1,NW,NG)+DTAUKV(K,NW,NG)
353        END DO
354     END DO           
355  END DO                 ! end full gauss loop
356
357
358!  Titan's outputs (J.V.O, 2016)===============================================
359!      do l=1,L_NLAYRAD
360!         do nw=1,L_NSPECTV
361!          INT_DTAU(L,NW) = 0.0d+0
362!            DO NG=1,L_NGAUSS
363!               INT_DTAU(L,NW)= INT_DTAU(L,NW) + dtauv(L,nw,ng)*gweight(NG)
364!            enddo
365!         enddo
366!      enddo
367
368!       do nw=1,L_NSPECTV
369!          write(str2,'(i2.2)') nw
370!         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))
371!          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))       
372!       enddo 
373
374! ============================================================================== 
375
376
377  return
378
379
380end subroutine optcv
Note: See TracBrowser for help on using the repository browser.