Ignore:
Timestamp:
Jun 8, 2017, 1:03:06 AM (7 years ago)
Author:
jvatant
Message:

Cleaning of the radiative code + coherence between optci and optcv :

+ Get rid of the actually dummy and confusing extra level L_LEVELS+1 for aerosols
+ Harmonisation between optci and optcv
+ Get rid of tauref (called nowhere) and csqtly of ini_radcommon subroutine

JVO

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.GENERIC/libf/phystd/optci.F90

    r1397 r1715  
    3232
    3333  real*8 DTAUI(L_NLAYRAD,L_NSPECTI,L_NGAUSS)
    34   real*8 DTAUKI(L_LEVELS+1,L_NSPECTI,L_NGAUSS)
     34  real*8 DTAUKI(L_LEVELS,L_NSPECTI,L_NGAUSS)
    3535  real*8 TAUI(L_NLEVRAD,L_NSPECTI,L_NGAUSS)
    3636  real*8 TAUCUMI(L_LEVELS,L_NSPECTI,L_NGAUSS)
     
    4242
    4343  ! for aerosols
    44   real*8  QXIAER(L_LEVELS+1,L_NSPECTI,NAERKIND)
    45   real*8  QSIAER(L_LEVELS+1,L_NSPECTI,NAERKIND)
    46   real*8  GIAER(L_LEVELS+1,L_NSPECTI,NAERKIND)
    47   real*8  TAUAERO(L_LEVELS+1,NAERKIND)
    48   real*8  TAUAEROLK(L_LEVELS+1,L_NSPECTI,NAERKIND)
     44  real*8  QXIAER(L_LEVELS,L_NSPECTI,NAERKIND)
     45  real*8  QSIAER(L_LEVELS,L_NSPECTI,NAERKIND)
     46  real*8  GIAER(L_LEVELS,L_NSPECTI,NAERKIND)
     47  real*8  TAUAERO(L_LEVELS,NAERKIND)
     48  real*8  TAUAEROLK(L_LEVELS,L_NSPECTI,NAERKIND)
    4949  real*8  TAEROS(L_LEVELS,L_NSPECTI,NAERKIND)
    5050
     
    135135  end do                    ! levels
    136136
    137 
     137  ! Spectral dependance of aerosol absorption
    138138  do iaer=1,naerkind
    139139     DO NW=1,L_NSPECTI
     
    147147
    148148     do K=2,L_LEVELS
    149 
    150 ! continuum absorption
    151         DCONT = 0.0d0
     149     
     150        DAERO=SUM(TAEROS(K,NW,1:naerkind)) ! aerosol absorption
     151
     152        DCONT = 0.0d0 ! continuum absorption
    152153
    153154        if(continuum.and.(.not.graybody))then
     
    225226        endif
    226227
    227 ! aerosol absorption
    228         DAERO=SUM(TAEROS(K,NW,1:naerkind))
    229 
    230228        do ng=1,L_NGAUSS-1
    231229
     
    286284  end do
    287285
    288   DTAUKI(L_LEVELS+1,1:L_NSPECTI,1:L_NGAUSS)=0.d0
    289286  !=======================================================================
    290287  !     Now the full treatment for the layers, where besides the opacity
     
    293290  do iaer=1,naerkind
    294291    DO NW=1,L_NSPECTI
    295      DO K=2,L_LEVELS+1
    296            TAUAEROLK(K,NW,IAER) = TAUAERO(K,IAER)*QSIAER(K,NW,IAER)
     292     DO K=2,L_LEVELS
     293           TAUAEROLK(K,NW,IAER) = TAUAERO(K,IAER)*QSIAER(K,NW,IAER) ! effect of scattering albedo
    297294     ENDDO
    298295    ENDDO
     
    300297 
    301298  DO NW=1,L_NSPECTI
    302      DO L=1,L_NLAYRAD
     299     DO L=1,L_NLAYRAD-1
    303300        K              = 2*L+1
    304301        btemp(L,NW) = SUM(TAUAEROLK(K,NW,1:naerkind)) + SUM(TAUAEROLK(K+1,NW,1:naerkind))
    305302     END DO ! L vertical loop
     303     
     304     ! Last level
     305     L           = L_NLAYRAD
     306     K           = 2*L+1   
     307     btemp(L,NW) = SUM(TAUAEROLK(K,NW,1:naerkind))
     308     
    306309  END DO                    ! NW spectral loop
    307310 
     
    309312  DO NW=1,L_NSPECTI
    310313     NG = L_NGAUSS
    311      DO L=1,L_NLAYRAD
     314     DO L=1,L_NLAYRAD-1
    312315
    313316        K              = 2*L+1
     
    334337
    335338     END DO ! L vertical loop
     339     
     340     ! Last level
     341     
     342     L              = L_NLAYRAD
     343     K              = 2*L+1
     344     DTAUI(L,nw,ng) = DTAUKI(K,NW,NG) ! + 1.e-50
     345
     346     atemp = 0.
     347     if(DTAUI(L,NW,NG) .GT. 1.0D-9) then
     348        do iaer=1,naerkind
     349           atemp = atemp + GIAER(K,NW,IAER)   * TAUAEROLK(K,NW,IAER)
     350        end do
     351        WBARI(L,nw,ng) = btemp(L,nw)  / DTAUI(L,NW,NG)
     352     else
     353        WBARI(L,nw,ng) = 0.0D0
     354        DTAUI(L,NW,NG) = 1.0D-9
     355     endif
     356
     357     if(btemp(L,nw) .GT. 0.0d0) then
     358        cosbi(L,NW,NG) = atemp/btemp(L,nw)
     359     else
     360        cosbi(L,NW,NG) = 0.0D0
     361     end if
     362     
    336363
    337364     ! Now the other Gauss points, if needed.
     
    340367        IF(TAUGSURF(NW,NG) .gt. TLIMIT) THEN
    341368
    342            DO L=1,L_NLAYRAD
     369           DO L=1,L_NLAYRAD-1
    343370              K              = 2*L+1
    344371              DTAUI(L,nw,ng) = DTAUKI(K,NW,NG)+DTAUKI(K+1,NW,NG)! + 1.e-50
     
    355382              cosbi(L,NW,NG) = cosbi(L,NW,L_NGAUSS)
    356383           END DO ! L vertical loop
     384           
     385           ! Last level
     386           L              = L_NLAYRAD
     387           K              = 2*L+1
     388           DTAUI(L,nw,ng) = DTAUKI(K,NW,NG)! + 1.e-50
     389
     390           if(DTAUI(L,NW,NG) .GT. 1.0D-9) then
     391
     392              WBARI(L,nw,ng) = btemp(L,nw)  / DTAUI(L,NW,NG)
     393
     394           else
     395              WBARI(L,nw,ng) = 0.0D0
     396              DTAUI(L,NW,NG) = 1.0D-9
     397           endif
     398
     399           cosbi(L,NW,NG) = cosbi(L,NW,L_NGAUSS)
     400           
    357401        END IF
    358402
Note: See TracChangeset for help on using the changeset viewer.