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

Last change on this file since 3090 was 3090, checked in by slebonnois, 13 months ago

BdeBatz? : Cleans microphysics and makes few corrections for physics

  • Property svn:executable set to *
File size: 21.6 KB
Line 
1SUBROUTINE OPTCV(PQMO,NLAY,ZLEV,PLEV,TMID,PMID,  &
2     DTAUV,TAUV,TAUCUMV,WBARV,COSBV,TAURAY,TAUGSURF,SEASHAZEFACT,&
3     POPT_HAZE,POPT_CLOUDS,CDCOLUMN)
4
5  use radinc_h
6  use radcommon_h, only: gasv,gasv_recomb,tlimit,Cmk,gzlat_ig, &
7                         tgasref,pfgasref,wnov,scalep,indv
8  use gases_h
9  use datafile_mod, only: haze_opt_file
10  use comcstfi_mod, only: pi,r
11  use callkeys_mod, only: continuum,graybody,callgasvis,corrk_recombin,   &
12                          callclouds,callmufi,seashaze,uncoupl_optic_haze,&
13                          opt4clouds,FHVIS, Fssa
14  use tracer_h, only: nmicro,nice,ices_indx
15
16  implicit none
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  !     Modified
29  !     --------
30  !     J. Vatant d'Ollone (2016-17)
31  !         --> Clean and adaptation to Titan
32  !     B. de Batz de Trenquelléon (2022-2023)
33  !         --> Clean and correction to Titan
34  !         --> New optics added for Titan's clouds
35  !     
36  !==================================================================
37  !     
38  !     THIS SUBROUTINE SETS THE OPTICAL CONSTANTS IN THE VISUAL 
39  !     IT CALCULATES FOR EACH LAYER, FOR EACH SPECTRAL INTERVAL IN THE VISUAL
40  !     LAYER: WBAR, DTAU, COSBAR
41  !     LEVEL: TAU
42  !     
43  !     TAUV(L,NW,NG) is the cumulative optical depth at the top of radiation code
44  !     layer L. NW is spectral wavelength interval, ng the Gauss point index.
45  !     
46  !     TLEV(L) - Temperature at the layer boundary
47  !     PLEV(L) - Pressure at the layer boundary (i.e. level)
48  !     GASV(NT,NPS,NW,NG) - Visible k-coefficients
49  !     
50  !-------------------------------------------------------------------
51
52
53  !==========================================================
54  ! Input/Output
55  !==========================================================
56  REAL*8, INTENT(IN)  :: PQMO(nlay,nmicro)  ! Tracers for microphysics optics (X/m2).
57  INTEGER, INTENT(IN) :: NLAY               ! Number of pressure layers (for pqmo)
58  REAL*8, INTENT(IN)  :: ZLEV(NLAY+1)
59  REAL*8, INTENT(IN)  :: PLEV(L_LEVELS)
60  REAL*8, INTENT(IN)  :: TMID(L_LEVELS), PMID(L_LEVELS)
61  REAL*8, INTENT(IN)  :: TAURAY(L_NSPECTV)
62  REAL*8, INTENT(IN)  :: SEASHAZEFACT(L_LEVELS)
63  INTEGER, INTENT(IN) :: CDCOLUMN
64 
65  REAL*8, INTENT(OUT) :: DTAUV(L_NLAYRAD,L_NSPECTV,L_NGAUSS)
66  REAL*8, INTENT(OUT) :: TAUV(L_NLEVRAD,L_NSPECTV,L_NGAUSS)
67  REAL*8, INTENT(OUT) :: TAUCUMV(L_LEVELS,L_NSPECTV,L_NGAUSS)
68  REAL*8, INTENT(OUT) :: COSBV(L_NLAYRAD,L_NSPECTV,L_NGAUSS)
69  REAL*8, INTENT(OUT) :: WBARV(L_NLAYRAD,L_NSPECTV,L_NGAUSS)
70  REAL*8, INTENT(OUT) :: TAUGSURF(L_NSPECTV,L_NGAUSS-1)
71  REAL*8, INTENT(OUT) :: POPT_HAZE(L_LEVELS,L_NSPECTI,3)
72  REAL*8, INTENT(OUT) :: POPT_CLOUDS(L_LEVELS,L_NSPECTI,3)
73  ! ==========================================================
74 
75  real*8 DTAUKV(L_LEVELS,L_NSPECTV,L_NGAUSS)
76
77  ! Titan customisation
78  ! J. Vatant d'Ollone (2016)
79  real*8 DHAZE_T(L_LEVELS,L_NSPECTV)
80  real*8 DHAZES_T(L_LEVELS,L_NSPECTV)
81  real*8 SSA_T(L_LEVELS,L_NSPECTV)
82  real*8 ASF_T(L_LEVELS,L_NSPECTV)
83  ! ==========================
84
85  integer L, NW, NG, K, LK, IAER
86  integer MT(L_LEVELS), MP(L_LEVELS), NP(L_LEVELS)
87  real*8  ANS, TAUGAS
88  real*8  TRAY(L_LEVELS,L_NSPECTV)
89  real*8  DPR(L_LEVELS), U(L_LEVELS)
90  real*8  LCOEF(4), LKCOEF(L_LEVELS,4)
91
92  real*8 DCONT
93  real*8 DRAYAER
94  double precision wn_cont, p_cont, p_air, T_cont, dtemp, dtempc
95  double precision p_cross
96
97  real*8  KCOEF(4)
98 
99  ! temporary variable to reduce memory access time to gasv
100  real*8 tmpk(2,2)
101
102  ! temporary variables for multiple aerosol calculation
103  real*8 atemp(L_NLAYRAD,L_NSPECTV)
104  real*8 btemp(L_NLAYRAD,L_NSPECTV)
105  real*8 ctemp(L_NLAYRAD,L_NSPECTV)
106
107  ! variables for k in units m^-1
108  real*8 dz(L_LEVELS)
109
110  integer igas, jgas, ilay
111
112  integer interm
113
114  ! Variables for haze optics
115  character(len=200) file_path
116  logical file_ok
117  integer dumch
118  real*8  dumwvl
119
120  ! Variables for new optics
121  integer iq, iw, FTYPE, CTYPE
122  real*8 m0as,m0af,m0ccn,m3as,m3af,m3ccn,m3cld
123  real*8 dtauaer_s,dtauaer_f,dtau_ccn,dtau_cld
124  real*8,save :: rhoaer_s(L_NSPECTV),ssa_s(L_NSPECTV),asf_s(L_NSPECTV)
125  real*8,save :: rhoaer_f(L_NSPECTV),ssa_f(L_NSPECTV),asf_f(L_NSPECTV)
126  real*8,save :: ssa_ccn(L_NSPECTI),asf_ccn(L_NSPECTI)
127  real*8,save :: ssa_cld(L_NSPECTV),asf_cld(L_NSPECTV)
128!$OMP THREADPRIVATE(rhoaer_s,rhoaer_f,ssa_s,ssa_f,ssa_cld,asf_s,asf_f,asf_cld)
129 
130  logical,save :: firstcall=.true.
131!$OMP THREADPRIVATE(firstcall)
132
133
134  !! AS: to save time in computing continuum (see bilinearbig)
135  IF (.not.ALLOCATED(indv)) THEN
136      ALLOCATE(indv(L_NSPECTV,ngasmx,ngasmx))
137      indv = -9999 ! this initial value means "to be calculated"
138  ENDIF
139 
140  ! Some initialisation beacause there's a pb with disr_haze at the limits (nw=1)
141  ! I should check this - For now we set vars to zero : better than nans - JVO 2017
142  DHAZE_T(:,:) = 0.0
143  SSA_T(:,:)   = 0.0
144  ASF_T(:,:)   = 0.0
145 
146  ! Load tabulated haze optical properties if needed.
147  ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
148  IF (firstcall .AND. callmufi .AND. (.NOT. uncoupl_optic_haze)) THEN
149     OPEN(12,file=TRIM(haze_opt_file),form='formatted') ! The file has been inquired in physiq_mod firstcall
150     READ(12,*) ! dummy header
151     DO NW=1,L_NSPECTI
152       READ(12,*) ! there's IR 1st
153     ENDDO
154     DO NW=1,L_NSPECTV
155       READ(12,*) dumch, dumwvl, rhoaer_f(nw), ssa_f(nw), asf_f(nw), rhoaer_s(nw), ssa_s(nw), asf_s(nw)
156     ENDDO
157     CLOSE(12)
158  ENDIF
159
160   !=======================================================================
161   !     Determine the total gas opacity throughout the column, for each
162   !     spectral interval, NW, and each Gauss point, NG.
163   !     Calculate the continuum opacities, i.e., those that do not depend on
164   !     NG, the Gauss index.
165
166   taugsurf(:,:) = 0.0
167   dpr(:)        = 0.0
168   lkcoef(:,:)   = 0.0
169
170   do K=2,L_LEVELS
171      ilay = L_NLAYRAD+1 - k/2 ! int. arithmetic => gives the gcm layer index (reversed)
172      DPR(k) = PLEV(K)-PLEV(K-1)
173
174      ! if we have continuum opacities, we need dz
175      dz(k) = dpr(k)*R*TMID(K)/(gzlat_ig(ilay)*PMID(K))
176      U(k)  = Cmk(ilay)*DPR(k)     ! only Cmk line in optcv.F     
177
178      call tpindex(PMID(K),TMID(K),pfgasref,tgasref,LCOEF,MT(K),MP(K))
179
180      do LK=1,4
181         LKCOEF(K,LK) = LCOEF(LK)
182      end do
183   end do ! L_LEVELS
184
185   ! Rayleigh scattering
186   do NW=1,L_NSPECTV
187      TRAY(1:4,NW)   = 1.d-30
188      do K=5,L_LEVELS
189         TRAY(K,NW)   = TAURAY(NW) * DPR(K)
190      end do ! L_LEVELS
191   end do
192 
193   do NW=1,L_NSPECTV
194      ! We ignore K=1...
195      do K=2,L_LEVELS
196         ! int. arithmetic => gives the gcm layer index (reversed)
197         ilay = L_NLAYRAD+1 - k/2
198         
199         ! Optics coupled with the microphysics :
200         IF (callmufi .AND. (.NOT. uncoupl_optic_haze)) THEN
201
202            !==========================================================================================
203            ! Old optics (must have callclouds = .false.):
204            !==========================================================================================
205            IF (.NOT. opt4clouds) THEN
206               m3as = pqmo(ilay,2) / 2.0
207               m3af = pqmo(ilay,4) / 2.0
208               ! Cut-off (here for p = 2.7e3Pa / alt = 70km)
209               IF (ilay .lt. 23) THEN
210                  m3as = pqmo(23,2) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(24)-zlev(23))
211                  m3af = pqmo(23,4) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(24)-zlev(23))
212               ENDIF
213               
214               dtauaer_s = m3as*rhoaer_s(nw)
215               dtauaer_f = m3af*rhoaer_f(nw)
216
217            !==========================================================================================
218            ! New optics :
219            !==========================================================================================
220            ELSE
221               iw = (L_NSPECTV + 1) - NW ! Visible first and return
222               !-----------------------------
223               ! HAZE (Spherical + Fractal) :
224               !-----------------------------
225               FTYPE = 1
226
227               ! Spherical aerosols :
228               !---------------------
229               CTYPE = 5
230               m0as  = pqmo(ilay,1) / 2.0
231               m3as  = pqmo(ilay,2) / 2.0
232               ! If not callclouds : must have a cut-off (here for p = 2.7e3Pa / alt = 70km)
233               IF (.NOT. callclouds) THEN
234                  IF (ilay .lt. 23) THEN
235                     m0as = pqmo(23,1) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(24)-zlev(23))
236                     m3as = pqmo(23,2) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(24)-zlev(23))
237                  ENDIF
238               ENDIF
239               call get_haze_and_cloud_opacity(FTYPE,CTYPE,m0as,m3as,iw,dtauaer_s,ssa_s(nw),asf_s(nw))
240
241               ! Fractal aerosols :
242               !-------------------
243               CTYPE = FTYPE
244               m0af  = pqmo(ilay,3) / 2.0
245               m3af  = pqmo(ilay,4) / 2.0
246               ! If not callclouds : must have a cut-off (here for p = 2.7e3Pa / alt = 70km)
247               IF (.NOT. callclouds) THEN
248                  IF (ilay .lt. 23) THEN
249                     m0af = pqmo(23,3) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(24)-zlev(23))
250                     m3af = pqmo(23,4) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(24)-zlev(23))
251                  ENDIF
252               ENDIF
253               call get_haze_and_cloud_opacity(FTYPE,CTYPE,m0af,m3af,iw,dtauaer_f,ssa_f(nw),asf_f(nw))
254            ENDIF
255
256            ! Total of Haze opacity (dtau), SSA (w) and ASF (COS) :
257            DHAZE_T(k,nw) = dtauaer_s + dtauaer_f
258            IF (dtauaer_s + dtauaer_f .GT. 1.D-30) THEN
259               SSA_T(k,nw) = ( dtauaer_s*ssa_s(nw) + dtauaer_f*ssa_f(nw) ) / ( dtauaer_s+dtauaer_f )
260               ASF_T(k,nw) = ( dtauaer_s*ssa_s(nw)*asf_s(nw) + dtauaer_f*ssa_f(nw)*asf_f(nw) )  &
261                             / ( ssa_s(nw)*dtauaer_s + ssa_f(nw)*dtauaer_f )
262            ELSE
263               DHAZE_T(k,nw) = 0.D0
264               SSA_T(k,nw)   = 1.0
265               ASF_T(k,nw)   = 1.0
266            ENDIF
267            ! Diagnostics for the haze :
268            POPT_HAZE(k,nw,1) = DHAZE_T(k,nw) ! dtau
269            POPT_HAZE(k,nw,2) = SSA_T(k,nw)   ! wbar
270            POPT_HAZE(k,nw,3) = ASF_T(k,nw)   ! gbar
271
272            !---------------------
273            ! CLOUDS (Spherical) :
274            !---------------------
275            IF (callclouds) THEN
276               CTYPE = 0
277               m0ccn = pqmo(ilay,5) / 2.0
278               m3ccn = pqmo(ilay,6) / 2.0
279               m3cld = 0.0d0
280               
281               ! Clear / Dark column method :
282               !-----------------------------
283
284               ! CCN's SSA :
285               call get_haze_and_cloud_opacity(FTYPE,FTYPE,m0ccn,m3ccn,iw,dtau_ccn,ssa_ccn(nw),asf_ccn(nw))
286
287               ! Clear column (CCN, C2H2, C2H6, HCN) :
288               IF (CDCOLUMN == 0) THEN
289                  DO iq = 2, nice
290                     m3cld  = m3cld + (pqmo(ilay,ices_indx(iq)) / 2.0)
291                  ENDDO
292                  call get_haze_and_cloud_opacity(FTYPE,CTYPE,m0ccn,m3cld,iw,dtau_cld,ssa_cld(nw),asf_cld(nw))
293               
294               ! Dark column (CCN, CH4, C2H2, C2H6, HCN) :
295               ELSEIF (CDCOLUMN == 1) THEN
296                  DO iq = 1, nice
297                     m3cld  = m3cld + (pqmo(ilay,ices_indx(iq)) / 2.0)
298                  ENDDO
299                  call get_haze_and_cloud_opacity(FTYPE,CTYPE,m0ccn,m3cld,iw,dtau_cld,ssa_cld(nw),asf_cld(nw))
300               
301               ELSE
302                  WRITE(*,*) 'WARNING OPTCV.F90 : CDCOLUMN MUST BE 0 OR 1'
303                  WRITE(*,*) 'WE USE DARK COLUMN ...'
304                  DO iq = 1, nice
305                     m3cld  = m3cld + (pqmo(ilay,ices_indx(iq)) / 2.0)
306                  ENDDO
307                  call get_haze_and_cloud_opacity(FTYPE,CTYPE,m0ccn,m3cld,iw,dtau_cld,ssa_cld(nw),asf_cld(nw))
308               ENDIF
309
310               ! For small dropplets, opacity of nucleus dominates
311               ssa_cld(nw) = (ssa_ccn(nw)*m3ccn + ssa_cld(nw)*m3cld) / (m3ccn + m3cld)
312               ssa_cld(nw) = ssa_cld(nw) * Fssa
313
314               ! Total of Haze + Clouds opacity (dtau), SSA (w) and ASF (COS) :
315               DHAZE_T(k,nw) = dtauaer_s + dtauaer_f + dtau_cld
316               IF (DHAZE_T(k,nw) .GT. 1.D-30) THEN
317                  SSA_T(k,nw) = ( dtauaer_s*ssa_s(nw) + dtauaer_f*ssa_f(nw) + dtau_cld*ssa_cld(nw) ) / ( dtauaer_s+dtauaer_f+dtau_cld )
318                  ASF_T(k,nw) = ( dtauaer_s*ssa_s(nw)*asf_s(nw) + dtauaer_f*ssa_f(nw)*asf_f(nw) + dtau_cld*ssa_cld(nw)*asf_cld(nw) ) &
319                                / ( ssa_s(nw)*dtauaer_s + ssa_f(nw)*dtauaer_f + ssa_cld(nw)*dtau_cld )
320               ELSE
321                  DHAZE_T(k,nw) = 0.D0
322                  SSA_T(k,nw)   = 1.0
323                  ASF_T(k,nw)   = 1.0
324               ENDIF
325               
326               ! Tuning of optical properties (now useless...) :
327               DHAZE_T(k,nw) = DHAZE_T(k,nw) * (FHVIS * (1-SSA_T(k,nw)) + SSA_T(k,nw))
328               SSA_T(k,nw)   = SSA_T(k,nw) / (FHVIS * (1-SSA_T(k,nw)) + SSA_T(k,nw))
329
330               ! Diagnostics for clouds :
331               POPT_CLOUDS(k,nw,1) = DHAZE_T(k,nw) ! dtau
332               POPT_CLOUDS(k,nw,2) = SSA_T(k,nw)   ! wbar
333               POPT_CLOUDS(k,nw,3) = ASF_T(k,nw)   ! gbar
334           
335            ELSE
336               ! Diagnostics for clouds :
337               POPT_CLOUDS(k,nw,1) = 0.D0 ! dtau
338               POPT_CLOUDS(k,nw,2) = 1.0  ! wbar
339               POPT_CLOUDS(k,nw,3) = 1.0  ! gbar
340            ENDIF
341           
342         ! Optics and microphysics no coupled :
343         ELSE
344            ! Call fixed vertical haze profile of extinction - same for all columns
345            call disr_haze(dz(k),plev(k),wnov(nw),DHAZE_T(k,nw),SSA_T(k,nw),ASF_T(k,nw))
346            if (seashaze) DHAZE_T(k,nw) = DHAZE_T(k,nw)*seashazefact(k)
347            ! Diagnostics for the haze :
348            POPT_HAZE(k,nw,1) = DHAZE_T(k,nw) ! dtau
349            POPT_HAZE(k,nw,2) = SSA_T(k,nw)   ! wbar
350            POPT_HAZE(k,nw,3) = ASF_T(k,nw)   ! gbar
351            ! Diagnostics for clouds :
352            POPT_CLOUDS(k,nw,1) = 0.D0 ! dtau
353            POPT_CLOUDS(k,nw,2) = 1.0  ! wbar
354            POPT_CLOUDS(k,nw,3) = 1.0  ! gbar
355         ENDIF ! ENDIF callmufi
356         
357         !JL18 It seems to be good to have aerosols in the first "radiative layer" of the gcm in the IR
358         !   but visible does not handle very well diffusion in first layer.
359         !   The tauaero and tauray are thus set to 0 (a small value for rayleigh because the code crashes otherwise)
360         !   in the 4 first semilayers in optcv, but not optci.
361         !   This solves random variations of the sw heating at the model top.
362         if (k<5)  DHAZE_T(K,:) = 0.0
363         
364         DRAYAER = TRAY(K,NW)
365         ! DRAYAER is Tau RAYleigh scattering, plus AERosol opacity
366         DRAYAER = DRAYAER + DHAZE_T(K,NW) ! Titan's aerosol
367
368         DCONT = 0.0 ! continuum absorption
369
370         if(continuum.and.(.not.graybody).and.callgasvis)then
371            ! include continua if necessary
372            wn_cont = dble(wnov(nw))
373            T_cont  = dble(TMID(k))
374            do igas=1,ngasmx
375
376               p_cont  = dble(PMID(k)*scalep*gfrac(igas,ilay))
377
378               dtemp=0.0
379               if(igas.eq.igas_N2)then
380
381                  interm = indv(nw,igas,igas)
382   !                 call interpolateN2N2(wn_cont,T_cont,p_cont,dtemp,.false.,interm)
383                  indv(nw,igas,igas) = interm
384                  ! only goes to 500 cm^-1, so unless we're around a cold brown dwarf, this is irrelevant in the visible
385
386               elseif(igas.eq.igas_H2)then
387
388                  ! first do self-induced absorption
389                  interm = indv(nw,igas,igas)
390                  call interpolateH2H2(wn_cont,T_cont,p_cont,dtemp,.false.,interm)
391                  indv(nw,igas,igas) = interm
392
393                  ! then cross-interactions with other gases
394                  do jgas=1,ngasmx
395                     p_cross = dble(PMID(k)*scalep*gfrac(jgas,ilay))
396                     dtempc  = 0.0
397                     if(jgas.eq.igas_N2)then
398                        interm = indv(nw,igas,jgas)
399                        call interpolateN2H2(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm)
400                        indv(nw,igas,jgas) = interm
401                        ! should be irrelevant in the visible
402                     endif
403                     dtemp = dtemp + dtempc
404                  enddo
405
406                  elseif(igas.eq.igas_CH4)then
407
408                  ! first do self-induced absorption
409                  interm = indv(nw,igas,igas)
410                  call interpolateCH4CH4(wn_cont,T_cont,p_cont,dtemp,.false.,interm)
411                  indv(nw,igas,igas) = interm
412
413                  ! then cross-interactions with other gases
414                  do jgas=1,ngasmx
415                     p_cross = dble(PMID(k)*scalep*gfrac(jgas,ilay))
416                     dtempc  = 0.0
417                     if(jgas.eq.igas_N2)then
418                        interm = indv(nw,igas,jgas)
419                        call interpolateN2CH4(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm)
420                        indv(nw,igas,jgas) = interm
421                     endif
422                     dtemp = dtemp + dtempc
423                  enddo
424
425               endif
426
427               DCONT = DCONT + dtemp
428
429            enddo
430
431            DCONT = DCONT*dz(k)
432
433         endif
434
435         do ng=1,L_NGAUSS-1
436
437            ! Now compute TAUGAS
438
439            ! JVO 2017 : added tmpk because the repeated calls to gasi/v increased dramatically
440            ! the execution time of optci/v -> ~ factor 2 on the whole radiative
441            ! transfer on the tested simulations !
442
443            if (corrk_recombin) then
444               tmpk = GASV_RECOMB(MT(K):MT(K)+1,MP(K):MP(K)+1,NW,NG)
445            else
446               tmpk = GASV(MT(K):MT(K)+1,MP(K):MP(K)+1,1,NW,NG)
447            endif
448               
449            KCOEF(1) = tmpk(1,1) ! KCOEF(1) = GASV(MT(K),MP(K),1,NW,NG)
450            KCOEF(2) = tmpk(1,2) ! KCOEF(2) = GASV(MT(K),MP(K)+1,1,NW,NG)
451            KCOEF(3) = tmpk(2,2) ! KCOEF(3) = GASV(MT(K)+1,MP(K)+1,1,NW,NG)
452            KCOEF(4) = tmpk(2,1) ! KCOEF(4) = GASV(MT(K)+1,MP(K),1,NW,NG)
453
454            ! Interpolate the gaseous k-coefficients to the requested T,P values
455
456            ANS = LKCOEF(K,1)*KCOEF(1) + LKCOEF(K,2)*KCOEF(2) +            &
457                  LKCOEF(K,3)*KCOEF(3) + LKCOEF(K,4)*KCOEF(4)
458
459
460            TAUGAS  = U(k)*ANS
461
462            TAUGSURF(NW,NG) = TAUGSURF(NW,NG) + TAUGAS + DCONT
463            DTAUKV(K,nw,ng) = TAUGAS &
464                              + DRAYAER & ! DRAYAER includes all scattering contributions
465                              + DCONT ! For parameterized continuum aborption
466         end do
467
468         ! Now fill in the "clear" part of the spectrum (NG = L_NGAUSS),
469         ! which holds continuum opacity only
470
471         NG              = L_NGAUSS
472         DTAUKV(K,nw,ng) = DRAYAER + DCONT ! Scattering + parameterized continuum absorption, including Titan's haze
473
474      end do ! k = L_LEVELS
475   end do ! nw = L_NSPECTI
476
477  !=======================================================================
478  !     Now the full treatment for the layers, where besides the opacity
479  !     we need to calculate the scattering albedo and asymmetry factors
480  ! ======================================================================
481
482  ! Haze scattering
483            !JL18 It seems to be good to have aerosols in the first "radiative layer" of the gcm in the IR
484            !   but not in the visible
485            !   The dhaze_s is thus set to 0 in the 4 first semilayers in optcv, but not optci.
486            !   This solves random variations of the sw heating at the model top.
487  DO NW=1,L_NSPECTV
488    DHAZES_T(1:4,NW) = 0.d0
489    DO K=5,L_LEVELS
490      DHAZES_T(K,NW) = DHAZE_T(K,NW) * SSA_T(K,NW) ! effect of scattering albedo on haze
491    ENDDO
492  ENDDO
493
494  ! NW spectral loop
495  DO NW=1,L_NSPECTV
496   ! L vertical loop
497   DO L=1,L_NLAYRAD-1
498      K = 2*L+1
499           atemp(L,NW) = ASF_T(K,NW)*DHAZES_T(K,NW) + ASF_T(K+1,NW)*DHAZES_T(K+1,NW)
500      btemp(L,NW) = DHAZES_T(K,NW) + DHAZES_T(K+1,NW)
501           ctemp(L,NW) = btemp(L,NW) + 0.9999*(TRAY(K,NW) + TRAY(K+1,NW)) ! JVO 2017 : does this 0.999 is really meaningful ?
502           btemp(L,NW) = btemp(L,NW) + TRAY(K,NW) + TRAY(K+1,NW)
503      COSBV(L,NW,1:L_NGAUSS) = atemp(L,NW)/btemp(L,NW)
504   END DO
505     
506   ! Last level
507   L = L_NLAYRAD
508   K = 2*L+1
509   atemp(L,NW) = ASF_T(K,NW)*DHAZES_T(K,NW)
510   btemp(L,NW) = DHAZES_T(K,NW)
511   ctemp(L,NW) = btemp(L,NW) + 0.9999*TRAY(K,NW) ! JVO 2017 : does this 0.999 is really meaningful ?
512   btemp(L,NW) = btemp(L,NW) + TRAY(K,NW)
513   COSBV(L,NW,1:L_NGAUSS) = atemp(L,NW)/btemp(L,NW)
514  END DO
515
516  ! NG Gauss loop
517  DO NG=1,L_NGAUSS
518   ! NW spectral loop
519   DO NW=1,L_NSPECTV
520      ! L vertical loop
521      DO L=1,L_NLAYRAD-1
522         K = 2*L+1
523         DTAUV(L,nw,ng) = DTAUKV(K,NW,NG) + DTAUKV(K+1,NW,NG)
524         WBARV(L,nw,ng) = ctemp(L,NW) / DTAUV(L,nw,ng)
525      END DO
526
527      ! Last level
528      L = L_NLAYRAD
529      K = 2*L+1
530           DTAUV(L,nw,ng) = DTAUKV(K,NW,NG)
531      WBARV(L,NW,NG) = ctemp(L,NW) / DTAUV(L,NW,NG)
532   END DO
533  END DO
534
535  ! Total extinction optical depths
536  !DO NG=1,L_NGAUSS       ! full gauss loop
537  !   DO NW=1,L_NSPECTV       
538  !      TAUCUMV(1,NW,NG)=0.0D0
539  !      DO K=2,L_LEVELS
540  !         TAUCUMV(K,NW,NG)=TAUCUMV(K-1,NW,NG)+DTAUKV(K,NW,NG)
541  !      END DO
542  !      DO L=1,L_NLAYRAD
543  !         TAUV(L,NW,NG)=TAUCUMV(2*L,NW,NG)
544  !      END DO
545  !      TAUV(L,NW,NG)=TAUCUMV(2*L_NLAYRAD+1,NW,NG)
546  !   END DO           
547  !END DO                 ! end full gauss loop
548 
549  TAUCUMV(:,:,:) = DTAUKV(:,:,:)
550  DO L=1,L_NLAYRAD
551   TAUV(L,:,:)=TAUCUMV(2*L,:,:)
552  END DO
553  TAUV(L,:,:)=TAUCUMV(2*L_NLAYRAD+1,:,:)
554
555  if(firstcall) firstcall = .false.
556
557  return
558
559
560end subroutine optcv
Note: See TracBrowser for help on using the repository browser.