Ignore:
Timestamp:
Aug 20, 2025, 4:25:12 PM (4 days ago)
Author:
emillour
Message:

Mars PCM:
Some code tidying:

  • turn aeroptproperties, albedocaps, cvmgp and convadj into modules
  • remove useless check in callradite
  • clean module vlz_fi (remove obsolete #ifdef CRAY cpp alternatives)
  • remove function cvmgp since it was only called under #ifdef CRAY in vlz_fi

EM

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/libf/phymars/swr_toon.F

    r3740 r3901  
    1515      use callkeys_mod, only: rayleigh
    1616      use swrayleigh_mod, only: swrayleigh
     17      use cvmgt_mod, only: cvmgt
    1718     
    1819      IMPLICIT NONE
     
    7576C     ARGUMENTS
    7677C     ---------
    77       INTEGER KDLON, KFLEV, KNU
    78       REAL aerosol(NDLO2,KFLEV,naerkind), albedo(NDLO2,2),
    79      S    PDSIG(NDLO2,KFLEV),PSEC(NDLO2)
    80 
    81       REAL QVISsQREF3d(NDLO2,KFLEV,nsun,naerkind) 
    82       REAL omegaVIS3d(NDLO2,KFLEV,nsun,naerkind)   
    83       REAL gVIS3d(NDLO2,KFLEV,nsun,naerkind)
    84 
    85       REAL PPSOL(NDLO2)
    86       REAL PFD(NDLO2,KFLEV+1),PFU(NDLO2,KFLEV+1)
    87       REAL PRMU(NDLO2)
     78      INTEGER,INTENT(IN) :: KDLON, KFLEV, KNU
     79      REAL,INTENT(IN) :: aerosol(NDLO2,KFLEV,naerkind), albedo(NDLO2,2)
     80      REAL,INTENT(IN) :: PDSIG(NDLO2,KFLEV),PSEC(NDLO2)
     81
     82      REAL,INTENT(IN) :: QVISsQREF3d(NDLO2,KFLEV,nsun,naerkind) 
     83      REAL,INTENT(IN) :: omegaVIS3d(NDLO2,KFLEV,nsun,naerkind)   
     84      REAL,INTENT(IN) :: gVIS3d(NDLO2,KFLEV,nsun,naerkind)
     85
     86      REAL,INTENT(IN) :: PPSOL(NDLO2)
     87      REAL,INTENT(OUT) :: PFD(NDLO2,KFLEV+1),PFU(NDLO2,KFLEV+1)
     88      REAL,INTENT(IN) :: PRMU(NDLO2)
    8889
    8990C     LOCAL ARRAYS
     
    104105c     End part added by Tran The Trung
    105106
    106 c     Function
    107 c     --------
    108       real CVMGT
    109 
    110107c Computing TOTAL single scattering parameters by adding
    111108c  properties of all the NAERKIND kind of aerosols
     
    117114            ZTAUAZ(JL,JK) = 0.
    118115         END DO
    119          DO 106 JAE=1,naerkind
    120             DO 105 JL = 1 , KDLON
     116         DO JAE=1,naerkind
     117            DO JL = 1 , KDLON
    121118c              Mean Extinction optical depth in the spectral band
    122119c              ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
     
    133130     S           QVISsQREF3d(JL,JK,KNU,JAE)*
    134131     &           omegaVIS3d(JL,JK,KNU,JAE)*gVIS3d(JL,JK,KNU,JAE)
    135  105        CONTINUE
    136  106     CONTINUE
     132            ENDDO
     133         ENDDO
    137134      END DO
    138135C     
    139136      DO JK = 1 , nlaylte
    140137         DO JL = 1 , KDLON
     138           ! NB: function CVMGT(x1,x2,l) returns x1 if l==.true.
     139           !     or x2 otherwise. Maybe worth inlining someday?
    141140            ZCGAZ(JL,JK) = CVMGT( 0., ZCGAZ(JL,JK) / ZPIZAZ(JL,JK),
    142141     S            (ZPIZAZ(JL,JK).EQ.0) )
Note: See TracChangeset for help on using the changeset viewer.