SUBROUTINE OPTCV(PQMO,NLAY,ZLEV,PLEV,TMID,PMID, & DTAUV,TAUV,TAUCUMV,WBARV,COSBV,TAURAY,TAUGSURF,SEASHAZEFACT,& DIAG_OPTH,DIAG_OPTT,CDCOLUMN) use radinc_h use radcommon_h, only: gasv,gasv_recomb,tlimit,Cmk,gzlat_ig, & tgasref,pfgasref,wnov,scalep,indv use gases_h use datafile_mod, only: haze_opt_file use comcstfi_mod, only: pi,r use callkeys_mod, only: continuum,graybody,callgasvis,corrk_recombin, & callclouds,callmufi,seashaze,uncoupl_optic_haze,& opt4clouds,FHVIS,FCVIS use tracer_h, only: nmicro,nice,ices_indx implicit none !================================================================== ! ! Purpose ! ------- ! Calculates shortwave optical constants at each level. ! ! Authors ! ------- ! Adapted from the NASA Ames code by R. Wordsworth (2009) ! ! Modified ! -------- ! J. Vatant d'Ollone (2016-17) ! --> Clean and adaptation to Titan ! B. de Batz de Trenquelléon (2022-2023) ! --> Clean and correction ! --> New optics added for Titan's clouds ! !================================================================== ! ! THIS SUBROUTINE SETS THE OPTICAL CONSTANTS IN THE VISUAL ! IT CALCULATES FOR EACH LAYER, FOR EACH SPECTRAL INTERVAL IN THE VISUAL ! LAYER: WBAR, DTAU, COSBAR ! LEVEL: TAU ! ! TAUV(L,NW,NG) is the cumulative optical depth at the top of radiation code ! layer L. NW is spectral wavelength interval, ng the Gauss point index. ! ! TLEV(L) - Temperature at the layer boundary ! PLEV(L) - Pressure at the layer boundary (i.e. level) ! GASV(NT,NPS,NW,NG) - Visible k-coefficients ! !------------------------------------------------------------------- !========================================================== ! Input/Output !========================================================== REAL*8, INTENT(IN) :: PQMO(nlay,nmicro) ! Tracers for microphysics optics (X/m2). INTEGER, INTENT(IN) :: NLAY ! Number of pressure layers (for pqmo) REAL*8, INTENT(IN) :: ZLEV(NLAY+1) REAL*8, INTENT(IN) :: PLEV(L_LEVELS) REAL*8, INTENT(IN) :: TMID(L_LEVELS), PMID(L_LEVELS) REAL*8, INTENT(IN) :: TAURAY(L_NSPECTV) REAL*8, INTENT(IN) :: SEASHAZEFACT(L_LEVELS) INTEGER, INTENT(IN) :: CDCOLUMN REAL*8, INTENT(OUT) :: DTAUV(L_NLAYRAD,L_NSPECTV,L_NGAUSS) REAL*8, INTENT(OUT) :: TAUV(L_NLEVRAD,L_NSPECTV,L_NGAUSS) REAL*8, INTENT(OUT) :: TAUCUMV(L_LEVELS,L_NSPECTV,L_NGAUSS) REAL*8, INTENT(OUT) :: COSBV(L_NLAYRAD,L_NSPECTV,L_NGAUSS) REAL*8, INTENT(OUT) :: WBARV(L_NLAYRAD,L_NSPECTV,L_NGAUSS) REAL*8, INTENT(OUT) :: TAUGSURF(L_NSPECTV,L_NGAUSS-1) REAL*8, INTENT(OUT) :: DIAG_OPTH(L_LEVELS,L_NSPECTV,6) REAL*8, INTENT(OUT) :: DIAG_OPTT(L_LEVELS,L_NSPECTV,6) ! ========================================================== real*8 DTAUKV(L_LEVELS,L_NSPECTV,L_NGAUSS) ! Titan customisation ! J. Vatant d'Ollone (2016) real*8 DHAZE_T(L_LEVELS,L_NSPECTV) real*8 DHAZES_T(L_LEVELS,L_NSPECTV) real*8 SSA_T(L_LEVELS,L_NSPECTV) real*8 ASF_T(L_LEVELS,L_NSPECTV) ! ========================== integer L, NW, NG, K, LK, IAER integer MT(L_LEVELS), MP(L_LEVELS), NP(L_LEVELS) real*8 ANS, TAUGAS real*8 TRAY(L_LEVELS,L_NSPECTV) real*8 DPR(L_LEVELS), U(L_LEVELS) real*8 LCOEF(4), LKCOEF(L_LEVELS,4) real*8 DCONT real*8 DRAYAER double precision wn_cont, p_cont, p_air, T_cont, dtemp, dtempc double precision p_cross real*8 KCOEF(4) ! temporary variable to reduce memory access time to gasv real*8 tmpk(2,2) ! temporary variables for multiple aerosol calculation real*8 atemp(L_NLAYRAD,L_NSPECTV) real*8 btemp(L_NLAYRAD,L_NSPECTV) real*8 ctemp(L_NLAYRAD,L_NSPECTV) ! variables for k in units m^-1 real*8 dz(L_LEVELS) integer igas, jgas, ilay integer interm ! Variables for haze optics character(len=200) file_path logical file_ok integer dumch real*8 dumwvl integer ilev_cutoff ! pressure level index of cutoff real*8 corr_haze real*8 corr_alb ! Variables for new optics integer iq, iw, FTYPE, CTYPE real*8 m0as,m0af,m0ccn,m3as,m3af,m3ccn,m3cld real*8 dtauaer_s,dtauaer_f,dtau_ccn,dtau_cld real*8,save :: rhoaer_s(L_NSPECTV),ssa_s(L_NSPECTV),asf_s(L_NSPECTV) real*8,save :: rhoaer_f(L_NSPECTV),ssa_f(L_NSPECTV),asf_f(L_NSPECTV) real*8,save :: ssa_ccn(L_NSPECTV),asf_ccn(L_NSPECTV) real*8,save :: ssa_cld(L_NSPECTV),asf_cld(L_NSPECTV) !$OMP THREADPRIVATE(rhoaer_s,rhoaer_f,ssa_s,ssa_f,ssa_cld,asf_s,asf_f,asf_cld) logical,save :: firstcall=.true. !$OMP THREADPRIVATE(firstcall) !! AS: to save time in computing continuum (see bilinearbig) IF (.not.ALLOCATED(indv)) THEN ALLOCATE(indv(L_NSPECTV,ngasmx,ngasmx)) indv = -9999 ! this initial value means "to be calculated" ENDIF ! Some initialisation beacause there's a pb with disr_haze at the limits (nw=1) ! I should check this - For now we set vars to zero : better than nans - JVO 2017 DHAZE_T(:,:) = 0.0 SSA_T(:,:) = 0.0 ASF_T(:,:) = 0.0 ! Load tabulated haze optical properties if needed. ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ IF (firstcall .AND. callmufi .AND. (.NOT. uncoupl_optic_haze)) THEN OPEN(12,file=TRIM(haze_opt_file),form='formatted') ! The file has been inquired in physiq_mod firstcall READ(12,*) ! dummy header DO NW=1,L_NSPECTI READ(12,*) ! there's IR 1st ENDDO DO NW=1,L_NSPECTV READ(12,*) dumch, dumwvl, rhoaer_f(nw), ssa_f(nw), asf_f(nw), rhoaer_s(nw), ssa_s(nw), asf_s(nw) ENDDO CLOSE(12) ENDIF !======================================================================= ! Determine the total gas opacity throughout the column, for each ! spectral interval, NW, and each Gauss point, NG. ! Calculate the continuum opacities, i.e., those that do not depend on ! NG, the Gauss index. taugsurf(:,:) = 0.0 dpr(:) = 0.0 lkcoef(:,:) = 0.0 ! Level of cutoff !~~~~~~~~~~~~~~~~ ilev_cutoff = 26 do K=2,L_LEVELS ilay = L_NLAYRAD+1 - k/2 ! int. arithmetic => gives the gcm layer index (reversed) DPR(k) = PLEV(K)-PLEV(K-1) ! if we have continuum opacities, we need dz dz(k) = dpr(k)*R*TMID(K)/(gzlat_ig(ilay)*PMID(K)) U(k) = Cmk(ilay)*DPR(k) ! only Cmk line in optcv.F call tpindex(PMID(K),TMID(K),pfgasref,tgasref,LCOEF,MT(K),MP(K)) do LK=1,4 LKCOEF(K,LK) = LCOEF(LK) end do end do ! L_LEVELS ! Rayleigh scattering do NW=1,L_NSPECTV TRAY(1:4,NW) = 1.d-30 do K=5,L_LEVELS TRAY(K,NW) = TAURAY(NW) * DPR(K) end do ! L_LEVELS end do DIAG_OPTH(:,:,:) = 0.D0 DIAG_OPTT(:,:,:) = 0.D0 do NW=1,L_NSPECTV ! We ignore K=1... do K=2,L_LEVELS ! int. arithmetic => gives the gcm layer index (reversed) ilay = L_NLAYRAD+1 - k/2 ! Optics coupled with the microphysics : IF (callmufi .AND. (.NOT. uncoupl_optic_haze)) THEN !========================================================================================== ! Old optics (must have callclouds = .false.): !========================================================================================== IF (.NOT. opt4clouds) THEN m3as = pqmo(ilay,2) / 2.0 m3af = pqmo(ilay,4) / 2.0 ! Cut-off (here for p = 2.7e3Pa / alt = 70km) IF (ilay .lt. ilev_cutoff) THEN m3as = pqmo(ilev_cutoff,2) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(ilev_cutoff+1)-zlev(ilev_cutoff)) m3af = pqmo(ilev_cutoff,4) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(ilev_cutoff+1)-zlev(ilev_cutoff)) ENDIF dtauaer_s = m3as*rhoaer_s(nw) dtauaer_f = m3af*rhoaer_f(nw) !========================================================================================== ! New optics : !========================================================================================== ELSE iw = (L_NSPECTV + 1) - NW ! Visible first and return !----------------------------- ! HAZE (Spherical + Fractal) : !----------------------------- FTYPE = 1 ! Spherical aerosols : !--------------------- CTYPE = 5 m0as = pqmo(ilay,1) / 2.0 m3as = pqmo(ilay,2) / 2.0 ! If not callclouds : must have a cut-off (here for p = 2.7e3Pa / alt = 70km) IF (.NOT. callclouds) THEN IF (ilay .lt. ilev_cutoff) THEN m0as = pqmo(ilev_cutoff,1) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(ilev_cutoff+1)-zlev(ilev_cutoff)) m3as = pqmo(ilev_cutoff,2) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(ilev_cutoff+1)-zlev(ilev_cutoff)) ENDIF ENDIF call get_haze_and_cloud_opacity(FTYPE,CTYPE,m0as,m3as,iw,dtauaer_s,ssa_s(nw),asf_s(nw)) ! Fractal aerosols : !------------------- CTYPE = FTYPE m0af = pqmo(ilay,3) / 2.0 m3af = pqmo(ilay,4) / 2.0 ! If not callclouds : must have a cut-off (here for p = 2.7e3Pa / alt = 70km) IF (.NOT. callclouds) THEN IF (ilay .lt. ilev_cutoff) THEN m0af = pqmo(ilev_cutoff,3) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(ilev_cutoff+1)-zlev(ilev_cutoff)) m3af = pqmo(ilev_cutoff,4) / 2.0 * (zlev(ilay+1)-zlev(ilay)) / (zlev(ilev_cutoff+1)-zlev(ilev_cutoff)) ENDIF ENDIF call get_haze_and_cloud_opacity(FTYPE,CTYPE,m0af,m3af,iw,dtauaer_f,ssa_f(nw),asf_f(nw)) ENDIF ! Tuning of optical properties for haze : !dtauaer_s = dtauaer_s * (FHVIS * (1-ssa_s(nw)) + ssa_s(nw)) !ssa_s(nw) = ssa_s(nw) / (FHVIS * (1-ssa_s(nw)) + ssa_s(nw)) dtauaer_f = dtauaer_f * (FHVIS * (1-ssa_f(nw)) + ssa_f(nw)) ssa_f(nw) = ssa_f(nw) / (FHVIS * (1-ssa_f(nw)) + ssa_f(nw)) ! Total of Haze opacity (dtau), SSA (w) and ASF (COS) : DHAZE_T(k,nw) = dtauaer_s + dtauaer_f IF (dtauaer_s + dtauaer_f .GT. 1.D-30) THEN SSA_T(k,nw) = ( dtauaer_s*ssa_s(nw) + dtauaer_f*ssa_f(nw) ) / ( dtauaer_s+dtauaer_f ) ASF_T(k,nw) = ( dtauaer_s*ssa_s(nw)*asf_s(nw) + dtauaer_f*ssa_f(nw)*asf_f(nw) ) & / ( ssa_s(nw)*dtauaer_s + ssa_f(nw)*dtauaer_f ) ELSE DHAZE_T(k,nw) = 0.D0 SSA_T(k,nw) = 1.0 ASF_T(k,nw) = 1.0 ENDIF ! Opacity and albedo adjustement below cutoff !---------------------------------------------- corr_haze=0.85-0.15*TANH((PMID(k)*100.-2500.)/250.) corr_alb =0.85-0.15*TANH((PMID(k)*100.-4500.)/1000.) IF (.NOT. callclouds) THEN IF (ilay .lt. ilev_cutoff) THEN DHAZE_T(k,nw) = DHAZE_T(k,nw) * corr_haze SSA_T(k,nw) = SSA_T(k,nw) * corr_alb ENDIF ENDIF ! Diagnostics for the haze : DIAG_OPTH(k,nw,1) = DHAZE_T(k,nw) ! dtau DIAG_OPTH(k,nw,2) = SSA_T(k,nw) ! wbar DIAG_OPTH(k,nw,3) = ASF_T(k,nw) ! gbar !--------------------- ! CLOUDS (Spherical) : !--------------------- IF (callclouds) THEN CTYPE = 0 m0ccn = pqmo(ilay,5) / 2.0 m3ccn = pqmo(ilay,6) / 2.0 m3cld = pqmo(ilay,6) / 2.0 ! Clear / Dark column method : !----------------------------- ! CCN's SSA : call get_haze_and_cloud_opacity(FTYPE,FTYPE,m0ccn,m3ccn,iw,dtau_ccn,ssa_ccn(nw),asf_ccn(nw)) ! Clear column (CCN, C2H2, C2H6, HCN, AC6H6) : IF (CDCOLUMN == 0) THEN DO iq = 2, nice m3cld = m3cld + (pqmo(ilay,ices_indx(iq)) / 2.0) ENDDO call get_haze_and_cloud_opacity(FTYPE,CTYPE,m0ccn,m3cld,iw,dtau_cld,ssa_cld(nw),asf_cld(nw)) ! Dark column (CCN, CH4, C2H2, C2H6, HCN, AC6H6) : ELSEIF (CDCOLUMN == 1) THEN DO iq = 1, nice m3cld = m3cld + (pqmo(ilay,ices_indx(iq)) / 2.0) ENDDO call get_haze_and_cloud_opacity(FTYPE,CTYPE,m0ccn,m3cld,iw,dtau_cld,ssa_cld(nw),asf_cld(nw)) ELSE WRITE(*,*) 'WARNING OPTCV.F90 : CDCOLUMN MUST BE 0 OR 1' WRITE(*,*) 'WE USE DARK COLUMN ...' DO iq = 1, nice m3cld = m3cld + (pqmo(ilay,ices_indx(iq)) / 2.0) ENDDO call get_haze_and_cloud_opacity(FTYPE,CTYPE,m0ccn,m3cld,iw,dtau_cld,ssa_cld(nw),asf_cld(nw)) ENDIF ! For small dropplets, opacity of nucleus dominates dtau_cld = (dtau_cld*m3ccn + dtau_cld*m3cld) / (m3ccn + m3cld) ssa_cld(nw) = (ssa_ccn(nw)*m3ccn + ssa_cld(nw)*m3cld) / (m3ccn + m3cld) ! Tuning of optical properties for clouds : dtau_cld = dtau_cld * (FCVIS * (1-ssa_cld(nw)) + ssa_cld(nw)) ssa_cld(nw) = ssa_cld(nw) / (FCVIS * (1-ssa_cld(nw)) + ssa_cld(nw)) ! Total of Haze + Clouds opacity (dtau), SSA (w) and ASF (COS) : DHAZE_T(k,nw) = dtauaer_s + dtauaer_f + dtau_cld IF (DHAZE_T(k,nw) .GT. 1.D-30) THEN 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 ) 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) ) & / ( ssa_s(nw)*dtauaer_s + ssa_f(nw)*dtauaer_f + ssa_cld(nw)*dtau_cld ) ELSE DHAZE_T(k,nw) = 0.D0 SSA_T(k,nw) = 1.0 ASF_T(k,nw) = 1.0 ENDIF ! Diagnostics for clouds : DIAG_OPTT(k,nw,1) = DHAZE_T(k,nw) ! dtau DIAG_OPTT(k,nw,2) = SSA_T(k,nw) ! wbar DIAG_OPTT(k,nw,3) = ASF_T(k,nw) ! gbar ELSE ! Diagnostics for clouds : DIAG_OPTT(k,nw,1) = 0.D0 ! dtau DIAG_OPTT(k,nw,2) = 1.0 ! wbar DIAG_OPTT(k,nw,3) = 1.0 ! gbar ENDIF ! Optics and microphysics no coupled : ELSE ! Call fixed vertical haze profile of extinction - same for all columns call disr_haze(dz(k),plev(k),wnov(nw),DHAZE_T(k,nw),SSA_T(k,nw),ASF_T(k,nw)) if (seashaze) DHAZE_T(k,nw) = DHAZE_T(k,nw)*seashazefact(k) ! Diagnostics for the haze : DIAG_OPTH(k,nw,1) = DHAZE_T(k,nw) ! dtau DIAG_OPTH(k,nw,2) = SSA_T(k,nw) ! wbar DIAG_OPTH(k,nw,3) = ASF_T(k,nw) ! gbar ! Diagnostics for clouds : DIAG_OPTT(k,nw,1) = 0.D0 ! dtau DIAG_OPTT(k,nw,2) = 1.0 ! wbar DIAG_OPTT(k,nw,3) = 1.0 ! gbar ENDIF ! ENDIF callmufi !JL18 It seems to be good to have aerosols in the first "radiative layer" of the gcm in the IR ! but visible does not handle very well diffusion in first layer. ! The tauaero and tauray are thus set to 0 (a small value for rayleigh because the code crashes otherwise) ! in the 4 first semilayers in optcv, but not optci. ! This solves random variations of the sw heating at the model top. if (k<5) DHAZE_T(K,:) = 0.0 DRAYAER = TRAY(K,NW) ! DRAYAER is Tau RAYleigh scattering, plus AERosol opacity DRAYAER = DRAYAER + DHAZE_T(K,NW) ! Titan's aerosol DCONT = 0.0 ! continuum absorption if(continuum.and.(.not.graybody).and.callgasvis)then ! include continua if necessary wn_cont = dble(wnov(nw)) T_cont = dble(TMID(k)) do igas=1,ngasmx p_cont = dble(PMID(k)*scalep*gfrac(igas,ilay)) dtemp=0.0 if(igas.eq.igas_N2)then interm = indv(nw,igas,igas) ! call interpolateN2N2(wn_cont,T_cont,p_cont,dtemp,.false.,interm) indv(nw,igas,igas) = interm ! only goes to 500 cm^-1, so unless we're around a cold brown dwarf, this is irrelevant in the visible elseif(igas.eq.igas_H2)then ! first do self-induced absorption interm = indv(nw,igas,igas) call interpolateH2H2(wn_cont,T_cont,p_cont,dtemp,.false.,interm) indv(nw,igas,igas) = interm ! then cross-interactions with other gases do jgas=1,ngasmx p_cross = dble(PMID(k)*scalep*gfrac(jgas,ilay)) dtempc = 0.0 if(jgas.eq.igas_N2)then interm = indv(nw,igas,jgas) call interpolateN2H2(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm) indv(nw,igas,jgas) = interm ! should be irrelevant in the visible endif dtemp = dtemp + dtempc enddo elseif(igas.eq.igas_CH4)then ! first do self-induced absorption interm = indv(nw,igas,igas) call interpolateCH4CH4(wn_cont,T_cont,p_cont,dtemp,.false.,interm) indv(nw,igas,igas) = interm ! then cross-interactions with other gases do jgas=1,ngasmx p_cross = dble(PMID(k)*scalep*gfrac(jgas,ilay)) dtempc = 0.0 if(jgas.eq.igas_N2)then interm = indv(nw,igas,jgas) call interpolateN2CH4(wn_cont,T_cont,p_cross,p_cont,dtempc,.false.,interm) indv(nw,igas,jgas) = interm endif dtemp = dtemp + dtempc enddo endif DCONT = DCONT + dtemp enddo DCONT = DCONT*dz(k) endif do ng=1,L_NGAUSS-1 ! Now compute TAUGAS ! JVO 2017 : added tmpk because the repeated calls to gasi/v increased dramatically ! the execution time of optci/v -> ~ factor 2 on the whole radiative ! transfer on the tested simulations ! if (corrk_recombin) then tmpk = GASV_RECOMB(MT(K):MT(K)+1,MP(K):MP(K)+1,NW,NG) else tmpk = GASV(MT(K):MT(K)+1,MP(K):MP(K)+1,1,NW,NG) endif KCOEF(1) = tmpk(1,1) ! KCOEF(1) = GASV(MT(K),MP(K),1,NW,NG) KCOEF(2) = tmpk(1,2) ! KCOEF(2) = GASV(MT(K),MP(K)+1,1,NW,NG) KCOEF(3) = tmpk(2,2) ! KCOEF(3) = GASV(MT(K)+1,MP(K)+1,1,NW,NG) KCOEF(4) = tmpk(2,1) ! KCOEF(4) = GASV(MT(K)+1,MP(K),1,NW,NG) ! Interpolate the gaseous k-coefficients to the requested T,P values ANS = LKCOEF(K,1)*KCOEF(1) + LKCOEF(K,2)*KCOEF(2) + & LKCOEF(K,3)*KCOEF(3) + LKCOEF(K,4)*KCOEF(4) TAUGAS = U(k)*ANS TAUGSURF(NW,NG) = TAUGSURF(NW,NG) + TAUGAS + DCONT DTAUKV(K,nw,ng) = TAUGAS & + DRAYAER & ! DRAYAER includes all scattering contributions + DCONT ! For parameterized continuum aborption end do ! Now fill in the "clear" part of the spectrum (NG = L_NGAUSS), ! which holds continuum opacity only NG = L_NGAUSS DTAUKV(K,nw,ng) = DRAYAER + DCONT ! Scattering + parameterized continuum absorption, including Titan's haze DIAG_OPTH(K,nw,4) = DRAYAER DIAG_OPTH(K,nw,5) = TAUGAS DIAG_OPTH(K,nw,6) = DCONT DIAG_OPTT(K,nw,4) = DRAYAER DIAG_OPTT(K,nw,5) = TAUGAS DIAG_OPTT(K,nw,6) = DCONT end do ! K = L_LEVELS end do ! nw = L_NSPECTV !======================================================================= ! Now the full treatment for the layers, where besides the opacity ! we need to calculate the scattering albedo and asymmetry factors ! ====================================================================== ! Haze scattering !JL18 It seems to be good to have aerosols in the first "radiative layer" of the gcm in the IR ! but not in the visible ! The dhaze_s is thus set to 0 in the 4 first semilayers in optcv, but not optci. ! This solves random variations of the sw heating at the model top. DO NW=1,L_NSPECTV DHAZES_T(1:4,NW) = 0.d0 DO K=5,L_LEVELS DHAZES_T(K,NW) = DHAZE_T(K,NW) * SSA_T(K,NW) ! effect of scattering albedo on haze ENDDO ENDDO ! NW spectral loop DO NW=1,L_NSPECTV ! L vertical loop DO L=1,L_NLAYRAD-1 K = 2*L+1 atemp(L,NW) = ASF_T(K,NW)*DHAZES_T(K,NW) + ASF_T(K+1,NW)*DHAZES_T(K+1,NW) btemp(L,NW) = DHAZES_T(K,NW) + DHAZES_T(K+1,NW) ctemp(L,NW) = btemp(L,NW) + 0.9999*(TRAY(K,NW) + TRAY(K+1,NW)) ! JVO 2017 : does this 0.999 is really meaningful ? btemp(L,NW) = btemp(L,NW) + TRAY(K,NW) + TRAY(K+1,NW) COSBV(L,NW,1:L_NGAUSS) = atemp(L,NW)/btemp(L,NW) END DO ! Last level L = L_NLAYRAD K = 2*L+1 atemp(L,NW) = ASF_T(K,NW)*DHAZES_T(K,NW) btemp(L,NW) = DHAZES_T(K,NW) ctemp(L,NW) = btemp(L,NW) + 0.9999*TRAY(K,NW) ! JVO 2017 : does this 0.999 is really meaningful ? btemp(L,NW) = btemp(L,NW) + TRAY(K,NW) COSBV(L,NW,1:L_NGAUSS) = atemp(L,NW)/btemp(L,NW) END DO ! NG Gauss loop DO NG=1,L_NGAUSS ! NW spectral loop DO NW=1,L_NSPECTV ! L vertical loop DO L=1,L_NLAYRAD-1 K = 2*L+1 DTAUV(L,nw,ng) = DTAUKV(K,NW,NG) + DTAUKV(K+1,NW,NG) WBARV(L,nw,ng) = ctemp(L,NW) / DTAUV(L,nw,ng) END DO ! Last level L = L_NLAYRAD K = 2*L+1 DTAUV(L,nw,ng) = DTAUKV(K,NW,NG) WBARV(L,NW,NG) = ctemp(L,NW) / DTAUV(L,NW,NG) END DO END DO ! Total extinction optical depths !DO NG=1,L_NGAUSS ! full gauss loop ! DO NW=1,L_NSPECTV ! TAUCUMV(1,NW,NG)=0.0D0 ! DO K=2,L_LEVELS ! TAUCUMV(K,NW,NG)=TAUCUMV(K-1,NW,NG)+DTAUKV(K,NW,NG) ! END DO ! DO L=1,L_NLAYRAD ! TAUV(L,NW,NG)=TAUCUMV(2*L,NW,NG) ! END DO ! TAUV(L,NW,NG)=TAUCUMV(2*L_NLAYRAD+1,NW,NG) ! END DO !END DO ! end full gauss loop TAUCUMV(:,:,:) = DTAUKV(:,:,:) DO L=1,L_NLAYRAD TAUV(L,:,:)=TAUCUMV(2*L,:,:) END DO TAUV(L,:,:)=TAUCUMV(2*L_NLAYRAD+1,:,:) if(firstcall) firstcall = .false. return end subroutine optcv