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

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

Useless argument Gweight in rad. tr. routines ( present in radcommon.h )
-JVO

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