Changeset 2677 for LMDZ5/trunk


Ignore:
Timestamp:
Oct 21, 2016, 1:52:36 PM (8 years ago)
Author:
oboucher
Message:

Reorganization of ocean_albedo.F90 calculation to speed up calculations.
Should be byte comparable. And speed up the routine by about 25-30%.

File:
1 edited

Legend:

Unmodified
Added
Removed
  • LMDZ5/trunk/libf/phylmd/ocean_albedo.F90

    r2227 r2677  
    33!
    44
    5 !     #########
    6 
    7 subroutine ocean_albedo(knon,zrmu0,knindex,pwind,SFRWL,alb_dir_new,alb_dif_new)
    8 
    9 
    10 !     ##################################################################
    11 !
     5SUBROUTINE ocean_albedo(knon,zrmu0,knindex,pwind,SFRWL,alb_dir_new,alb_dif_new)
     6!!
    127!!****  *ALBEDO_RS14* 
    138!!
    149!!    PURPOSE
    1510!!    -------
    16 !       computes the direct & diffuse albedo over open water
    17 !
    18 !     
     11!!     computes the direct & diffuse albedo over open water
     12!!     
    1913!!**  METHOD
    2014!!    ------
    21 !
     15!!
    2216!!    EXTERNAL
    2317!!    --------
     
    2519!!    IMPLICIT ARGUMENTS
    2620!!    ------------------
    27 !!
    2821!!     
    2922!!    REFERENCE
    3023!!    ---------
    31 !!
    3224!!     
    3325!!    AUTHOR
     
    3830!!    -------------
    3931!!      Original    03/2014
    40 !                   05/2014 R. Séférian & B. Decharme :: Adaptation to spectral
    41 !                   computation for diffuse and direct albedo
    42 !                   08/2014 S. Baek :: for wider wavelength range 200-4000nm and
    43 !                   adaptation to LMDZ + whitecap effect by Koepke + chrolophyll
    44 !                   map from climatology file
    45 !       
     32!!                  05/2014 R. Séférian & B. Decharme :: Adaptation to spectral
     33!!                  computation for diffuse and direct albedo
     34!!                  08/2014 S. Baek :: for wider wavelength range 200-4000nm and
     35!!                  adaptation to LMDZ + whitecap effect by Koepke + chrolophyll
     36!!                  map from climatology file
     37!!                  10/2016 O. Boucher :: some optimisation following R.
     38!!                  Seferian's work in the CNRM Model
     39!!       
    4640!-------------------------------------------------------------------------------
    4741!
     
    5044!
    5145USE ocean_albedo_para
    52 use dimphy
    53 !LF USE PARKIND1  ,ONLY : JPIM     ,JPRB
    54 use phys_state_var_mod, only : chl_con
     46USE dimphy
     47USE phys_state_var_mod, ONLY : chl_con
    5548!
    5649!
     
    6053!              -------------------------
    6154!
    62 
    6355include "clesphys.h"
    64 
    65 integer, intent(in) :: knon
    66 integer, dimension(klon), intent(in) :: knindex
    67 real, dimension(klon), intent(in) :: zrmu0,pwind
    68 real, dimension(klon,nsw), intent(out) ::  alb_dir_new,alb_dif_new
    69 real, dimension(6),intent(in) :: SFRWL
    70 
    71 
    72 !=== LOCAL VARIABLES
    73 
    74 REAL, parameter :: XPI=4.*atan(1.)
    75 
     56!
     57INTEGER, INTENT(IN) :: knon
     58INTEGER, DIMENSION(klon), INTENT(IN) :: knindex
     59REAL, DIMENSION(klon), INTENT(IN) :: zrmu0         !--cos(SZA) on full vector
     60REAL, DIMENSION(klon), INTENT(IN) :: pwind         !--wind speed on compressed vector
     61REAL, DIMENSION(6),INTENT(IN) :: SFRWL
     62REAL, DIMENSION(klon,nsw), INTENT(OUT) :: alb_dir_new, alb_dif_new
    7663!
    7764!*      0.2    declarations of local variables
    7865!              -------------------------
    7966!
    80 REAL, DIMENSION(klon)               :: ZCHL        ! surface chlorophyll
    81 REAL, DIMENSION(klon,NNWL)               :: ZDIR_ALB    ! direct  ocean surface albedo (spectral)
    82 REAL, DIMENSION(klon,NNWL)               :: ZSCA_ALB    ! diffuse ocean surface albedo (spectral)
    83 !
    84 INTEGER                         :: JI, JWL                  ! indexes
    85 REAL                            :: ZWL                      ! input parameter: wavelength and diffuse/direct fraction of light
    86 REAL:: ZSIG, ZREFM, ZXX2, ZR00, ZRR0, ZRRR                        ! computation variables
    87 REAL:: ZR22, ZUE, ZUE2, ZR11DF, ZALBT, ZFWC                       ! computation variables
    88 REAL:: ZCHLABS, ZAW, ZBW, ZAP, ZYLMD, ZBP550                      ! computation variables
    89 REAL:: ZBBP, ZNU, ZHB                                       ! computation variables
    90 REAL:: ZCOSZEN                                              ! Cosine of the zenith solar angle
    91 REAL:: ZR11, ZRW, ZRWDF, ZRDF                                   ! 4 components of the OSA
    92 ! new damping coefficient
    93 REAL:: ZDAMP     
    94 
     67REAL, DIMENSION(klon)           :: ZCHL        ! surface chlorophyll
     68REAL, DIMENSION(klon)           :: ZCOSZEN     ! Cosine of the zenith solar angle
     69!
     70INTEGER                         :: JWL, INU    ! indexes
     71REAL                            :: ZWL         ! input parameter: wavelength and diffuse/direct fraction of light
     72REAL:: ZCHLABS, ZAW, ZBW, ZREFM, ZYLMD         ! scalar computation variables
     73!
     74REAL, DIMENSION(klon) :: ZAP, ZXX2, ZR00, ZRR0, ZRRR               ! computation variables
     75REAL, DIMENSION(klon) :: ZR22, ZUE, ZUE2, ZR11DF, ZALBT            ! computation variables
     76REAL, DIMENSION(klon) :: ZBBP, ZNU, ZHB                            ! computation variables
     77REAL, DIMENSION(klon) :: ZR11, ZRW, ZRWDF, ZRDF                    ! 4 components of the OSA
     78REAL, DIMENSION(klon) :: ZSIG, ZFWC, ZWORK1, ZWORK2, ZWORK3
    9579!
    96 REAL            :: ZWORK                           ! dummy variable
    97 !
    98 !LF REAL(KIND=JPRB) :: ZHOOK_HANDLE
    99 !
    100 !-------------------------------------------------------------------------------
    101 !
    102 !
    103 
    104 
    105 
     80!--initialisations-------------
    10681!
    10782alb_dir_new(:,:) = 0.
    10883alb_dif_new(:,:) = 0.
    10984!
    110 ZDIR_ALB(:,:) = 0.
    111 ZSCA_ALB(:,:) = 0.
    112 !
    113 !
    114 
    115 !ZCHL(:) = CHL_CON!0.05 ! averaged global values for surface chlorophyll
    116 if(ok_chlorophyll)then
    117   do ji=1,knon
    118   ZCHL(ji)=CHL_CON(knindex(ji))
    119   enddo
    120 else
    121   ZCHL(:) = 0.05
    122 endif
    123 
    124 
    125 !
    126 DO JWL=1,NNWL           ! loop over the wavelength
    127 !
    128   DO JI=1,knon  ! loop over the grid points
    129 
    130 
    131     !---------------------------------------------------------------------------------
    132     ! 0- Compute baseline values
    133     !---------------------------------------------------------------------------------
     85! Initialisation of chlorophyll content
     86! ZCHL(:) = CHL_CON!0.05 ! averaged global values for surface chlorophyll
     87IF (ok_chlorophyll) THEN
     88  ZCHL(1:knon)=CHL_CON(knindex(1:knon))
     89ELSE
     90  ZCHL(1:knon) = 0.05
     91ENDIF
     92
     93! variables that do not depend on wavelengths
     94! loop over the grid points
     95! functions of chlorophyll content
     96ZWORK1(1:knon)= EXP(LOG(ZCHL(1:knon))*0.65)
     97ZWORK2(1:knon)= 0.416 * EXP(LOG(ZCHL(1:knon))*0.766)
     98ZWORK3(1:knon)= LOG10(ZCHL(1:knon))
     99! store the cosine of the solar zenith angle
     100ZCOSZEN(1:knon) = zrmu0(knindex(1:knon))
     101! Compute sigma derived from wind speed (Cox & Munk reflectance model)
     102ZSIG(1:knon)=SQRT(0.003+0.00512*PWIND(1:knon))
     103! original : correction for foam (Eq 16-17)
     104! has to be update once we have information from wave model (discussion with G. Madec)
     105ZFWC(1:knon)=3.97e-4*PWIND(1:knon)**1.59 ! Salisbury 2014 eq(2) at 37GHz, value in fraction
     106!
     107DO JWL=1,NNWL           ! loop over the wavelengths
     108  !
     109  !---------------------------------------------------------------------------------
     110  ! 0- Compute baseline values
     111  !---------------------------------------------------------------------------------
    134112   
    135     ! Get refractive index for the correspoding wavelength
    136     ZWL=XAKWL(JWL) !!!---------- wavelength value
    137     ZREFM= XAKREFM(JWL) !!!--------- refraction index value
     113  ! Get refractive index for the correspoding wavelength
     114  ZWL=XAKWL(JWL)      !!!--------- wavelength value
     115  ZREFM= XAKREFM(JWL) !!!--------- refraction index value
    138116 
    139 
    140     ! compute the cosine of the solar zenith angle
    141 !    ZCOSZEN = COS(XPI/2 - PZENITH(JI))
    142      ZCOSZEN = zrmu0(knindex(JI))
    143     ! Compute sigma derived from wind speed (Cox & Munk reflectance model)
    144     ZSIG=SQRT(0.003+0.00512*PWIND(JI))
    145 
    146    
    147     !---------------------------------------------------------------------------------
    148     ! 1- Compute direct surface albedo (ZR11)
    149     !---------------------------------------------------------------------------------
    150     !
    151     ZXX2=SQRT(1.0-(1.0-ZCOSZEN**2)/ZREFM**2)
    152     ZRR0=0.50*(((ZXX2-ZREFM*ZCOSZEN)/(ZXX2+ZREFM*ZCOSZEN))**2 +((ZCOSZEN-ZREFM*ZXX2)/(ZCOSZEN+ZREFM*ZXX2))**2)
    153     ZRRR=0.50*(((ZXX2-1.34*ZCOSZEN)/(ZXX2+1.34*ZCOSZEN))**2 +((ZCOSZEN-1.34*ZXX2)/(ZCOSZEN+1.34*ZXX2))**2)
    154     ZR11=ZRR0-(0.0152-1.7873*ZCOSZEN+6.8972*ZCOSZEN**2-8.5778*ZCOSZEN**3+4.071*ZSIG-7.6446*ZCOSZEN*ZSIG) &
    155         & * EXP(0.1643-7.8409*ZCOSZEN-3.5639*ZCOSZEN**2-2.3588*ZSIG+10.0538*ZCOSZEN*ZSIG)*ZRR0/ZRRR
    156     !
    157     !---------------------------------------------------------------------------------
    158     ! 2- Compute surface diffuse albedo (ZRDF)
    159     !---------------------------------------------------------------------------------
    160     ! Diffuse albedo from Jin et al., 2006 + estimation from diffuse fraction of
    161     ! light (relying later on AOD)
    162     ZRDF=-0.1482-0.012*ZSIG+0.1609*ZREFM-0.0244*ZSIG*ZREFM ! surface diffuse (Eq 5a-5b)
    163    
    164     !---------------------------------------------------------------------------------
    165     ! *- Determine absorption and backscattering
    166     ! coefficients to determine reflectance below the surface (Ro) once for all
    167     !
    168     ! *.1- Absorption by chlorophyll
    169     ZCHLABS= XAKACHL(JWL)
    170     ! *.2- Absorption by seawater
    171     ZAW= XAKAW3(JWL)
    172     ! *.3- Backscattering by seawater
    173     ZBW= XAKBW(JWL)
    174     ! *.4- Backscattering by chlorophyll
    175     ZYLMD = EXP(0.014*(440.0-ZWL))
    176     ZWORK= EXP(LOG(ZCHL(JI))*0.65)
    177     ZAP = 0.06*ZCHLABS*ZWORK +0.2*(XAW440+0.06*ZWORK)*ZYLMD
    178     ZBP550 = 0.416 * EXP(LOG(ZCHL(JI))*0.766)
    179    
    180     IF ( ZCHL(JI) > 2. ) THEN
    181       ZNU=0.
    182     ELSE
    183       IF ( ZCHL(JI) > 0.02 ) THEN
    184         ZWORK=LOG10(ZCHL(JI))
    185         ZNU=0.5*(ZWORK-0.3)
    186         ZBBP=(0.002+0.01*(0.5-0.25*ZWORK)*(ZWL/550.)**ZNU)*ZBP550
    187       ELSE
    188         ZBBP=0.019*(550./ZWL)*ZBP550       !ZBBPf=0.0113 at chl<=0.02
    189       ENDIF
     117  !---------------------------------------------------------------------------------
     118  ! 1- Compute direct surface albedo (ZR11)
     119  !---------------------------------------------------------------------------------
     120  !
     121  ZXX2(1:knon)=SQRT(1.0-(1.0-ZCOSZEN(1:knon)**2)/ZREFM**2)
     122  ZRR0(1:knon)=0.50*(((ZXX2(1:knon)-ZREFM*ZCOSZEN(1:knon))/(ZXX2(1:knon)+ZREFM*ZCOSZEN(1:knon)))**2 +  &
     123               ((ZCOSZEN(1:knon)-ZREFM*ZXX2(1:knon))/(ZCOSZEN(1:knon)+ZREFM*ZXX2(1:knon)))**2)
     124  ZRRR(1:knon)=0.50*(((ZXX2(1:knon)-1.34*ZCOSZEN(1:knon))/(ZXX2(1:knon)+1.34*ZCOSZEN(1:knon)))**2 + &
     125               ((ZCOSZEN(1:knon)-1.34*ZXX2(1:knon))/(ZCOSZEN(1:knon)+1.34*ZXX2(1:knon)))**2)
     126  ZR11(1:knon)=ZRR0(1:knon)-(0.0152-1.7873*ZCOSZEN(1:knon)+6.8972*ZCOSZEN(1:knon)**2-8.5778*ZCOSZEN(1:knon)**3+ &
     127               4.071*ZSIG(1:knon)-7.6446*ZCOSZEN(1:knon)*ZSIG(1:knon)) *  &
     128               EXP(0.1643-7.8409*ZCOSZEN(1:knon)-3.5639*ZCOSZEN(1:knon)**2-2.3588*ZSIG(1:knon)+ &
     129               10.0538*ZCOSZEN(1:knon)*ZSIG(1:knon))*ZRR0(1:knon)/ZRRR(1:knon)
     130  !
     131  !---------------------------------------------------------------------------------
     132  ! 2- Compute surface diffuse albedo (ZRDF)
     133  !---------------------------------------------------------------------------------
     134  ! Diffuse albedo from Jin et al., 2006 + estimation from diffuse fraction of
     135  ! light (relying later on AOD). CNRM model has opted for Eq 5b
     136  ZRDF(1:knon)=-0.1482-0.012*ZSIG(1:knon)+0.1609*ZREFM-0.0244*ZSIG(1:knon)*ZREFM ! surface diffuse (Eq 5a)
     137  !!ZRDF(1:knon)=-0.1479+0.1502*ZREFM-0.0176*ZSIG(1:knon)*ZREFM   ! surface diffuse (Eq 5b)
     138 
     139  !---------------------------------------------------------------------------------
     140  ! *- Determine absorption and backscattering
     141  ! coefficients to determine reflectance below the surface (Ro) once for all
     142  !
     143  ! *.1- Absorption by chlorophyll
     144  ZCHLABS= XAKACHL(JWL)
     145  ! *.2- Absorption by seawater
     146  ZAW= XAKAW3(JWL)
     147  ! *.3- Backscattering by seawater
     148  ZBW= XAKBW(JWL)
     149  ! *.4- Backscattering by chlorophyll
     150  ZYLMD = EXP(0.014*(440.0-ZWL))
     151  ZAP(1:knon) = 0.06*ZCHLABS*ZWORK1(1:knon) +0.2*(XAW440+0.06*ZWORK1(1:knon))*ZYLMD
     152   
     153  WHERE ( ZCHL(1:knon) > 0.02 )
     154    ZNU(:)=MIN(0.0,0.5*(ZWORK3(:)-0.3))
     155    ZBBP(:)=(0.002+0.01*(0.5-0.25*ZWORK3(:))*(ZWL/550.)**ZNU)*ZWORK2(:)
     156  ELSEWHERE
     157    ZBBP(:)=0.019*(550./ZWL)*ZWORK2(:)       !ZBBPf=0.0113 at chl<=0.02
     158  ENDWHERE
     159   
     160  ! Morel-Gentili(1991), Eq (12)
     161  ! ZHB=h/(h+2*ZBBPf*(1.-h))       
     162  ZHB(1:knon)=0.5*ZBW/(0.5*ZBW+ZBBP(1:knon))
     163   
     164  !---------------------------------------------------------------------------------
     165  ! 3- Compute direct water-leaving albedo (ZRW)
     166  !---------------------------------------------------------------------------------
     167  ! Based on Morel & Gentilli 1991 parametrization
     168  ZR22(1:knon)=0.48168549-0.014894708*ZSIG(1:knon)-0.20703885*ZSIG(1:knon)**2
     169
     170  ! Use Morel 91 formula to compute the direct reflectance
     171  ! below the surface
     172  ZR00(1:knon)=(0.5*ZBW+ZBBP(1:knon))/(ZAW+ZAP(1:knon)) *  &
     173               (0.6279-0.2227*ZHB(1:knon)-0.0513*ZHB(1:knon)**2 + &
     174               (-0.3119+0.2465*ZHB(1:knon))*ZCOSZEN(1:knon))
     175  ZRW(1:knon)=ZR00(1:knon)*(1.-ZR22(1:knon))/(1.-ZR00(1:knon)*ZR22(1:knon))
     176
     177  !---------------------------------------------------------------------------------
     178  ! 4- Compute diffuse water-leaving albedo (ZRWDF)
     179  !---------------------------------------------------------------------------------
     180  ! as previous water-leaving computation but assumes a uniform incidence of
     181  ! shortwave at surface (ue)
     182  ZUE=0.676               ! equivalent u_unif for diffuse incidence
     183  ZUE2=SQRT(1.0-(1.0-ZUE**2)/ZREFM**2)
     184  ZRR0(1:knon)=0.50*(((ZUE2-ZREFM*ZUE)/(ZUE2+ZREFM*ZUE))**2 +((ZUE-ZREFM*ZUE2)/(ZUE+ZREFM*ZUE2))**2)
     185  ZRRR(1:knon)=0.50*(((ZUE2-1.34*ZUE)/(ZUE2+1.34*ZUE))**2 +((ZUE-1.34*ZUE2)/(ZUE+1.34*ZUE2))**2)
     186  ZR11DF(1:knon)=ZRR0(1:knon)-(0.0152-1.7873*ZUE+6.8972*ZUE**2-8.5778*ZUE**3+4.071*ZSIG(1:knon)-7.6446*ZUE*ZSIG(1:knon)) * &
     187                 EXP(0.1643-7.8409*ZUE-3.5639*ZUE**2-2.3588*ZSIG(1:knon)+10.0538*ZUE*ZSIG(1:knon))*ZRR0(1:knon)/ZRRR(1:knon)
     188
     189  ! Use Morel 91 formula to compute the diffuse
     190  ! reflectance below the surface
     191  ZR00(1:knon)=(0.5*ZBW+ZBBP(1:knon))/(ZAW+ZAP(1:knon))*(0.6279-0.2227*ZHB(1:knon)-0.0513*ZHB(1:knon)**2 + (-0.3119+0.2465*ZHB(1:knon))*ZUE)
     192  ZRWDF(1:knon)=ZR00(1:knon)*(1.-ZR22(1:knon))*(1.-ZR11DF(1:knon))/(1.-ZR00(1:knon)*ZR22(1:knon))
     193   
     194  ! get waveband index inu for each nsw band
     195  SELECT CASE(nsw)
     196  CASE(2)
     197    IF (JWL.LE.49) THEN       ! from 200  to 680 nm
     198     inu=1
     199    ELSE                      ! from 690  to 4000 nm
     200     inu=2
    190201    ENDIF
    191    
    192     ! Morel-Gentili(1991), Eq (12)
    193     ! ZHB=h/(h+2*ZBBPf*(1.-h))       
    194     ZHB=0.5*ZBW/(0.5*ZBW+ZBBP)
    195    
    196     !---------------------------------------------------------------------------------
    197     ! 3- Compute direct water-leaving albedo (ZRW)
    198     !---------------------------------------------------------------------------------
    199     ! Based on Morel & Gentilli 1991 parametrization
    200     ZR22=0.48168549-0.014894708*ZSIG-0.20703885*ZSIG**2
    201     ! Use Morel 91 formula to compute the direct reflectance
    202     ! below the surface
    203     ZR00=(0.5*ZBW+ZBBP)/(ZAW+ZAP) *(0.6279-0.2227*ZHB-0.0513*ZHB**2 + (-0.3119+0.2465*ZHB)*ZCOSZEN)
    204     ZRW=ZR00*(1.-ZR22)*(1.-ZR11)/(1.-ZR00*ZR22)
    205 
    206     ZRW=ZR00*(1.-ZR22)/(1.-ZR00*ZR22)
    207     !---------------------------------------------------------------------------------
    208     ! 4- Compute diffuse water-leaving albedo (ZRWDF)
    209     !---------------------------------------------------------------------------------
    210     ! as previous water-leaving computation but assumes a uniform incidence of
    211     ! shortwave at surface (ue)
    212     ZUE=0.676               ! equivalent u_unif for diffuse incidence
    213     ZUE2=SQRT(1.0-(1.0-ZUE**2)/ZREFM**2)
    214     ZRR0=0.50*(((ZUE2-ZREFM*ZUE)/(ZUE2+ZREFM*ZUE))**2 +((ZUE-ZREFM*ZUE2)/(ZUE+ZREFM*ZUE2))**2)
    215     ZRRR=0.50*(((ZUE2-1.34*ZUE)/(ZUE2+1.34*ZUE))**2 +((ZUE-1.34*ZUE2)/(ZUE+1.34*ZUE2))**2)
    216     ZR11DF=ZRR0-(0.0152-1.7873*ZUE+6.8972*ZUE**2-8.5778*ZUE**3+4.071*ZSIG-7.6446*ZUE*ZSIG) &
    217           & * EXP(0.1643-7.8409*ZUE-3.5639*ZUE**2-2.3588*ZSIG+10.0538*ZUE*ZSIG)*ZRR0/ZRRR
    218     ! Use Morel 91 formula to compute the diffuse
    219     ! reflectance below the surface
    220     ZR00=(0.5*ZBW+ZBBP)/(ZAW+ZAP) *(0.6279-0.2227*ZHB-0.0513*ZHB**2 + (-0.3119+0.2465*ZHB)*ZUE)
    221     ZRWDF=ZR00*(1.-ZR22)*(1.-ZR11DF)/(1.-ZR00*ZR22)
    222    
    223     ! original : correction for foam (Eq 16-17)
    224     ZFWC=3.97e-4*PWIND(JI)**(1.59) ! Salisbury 2014 eq(2) at 37GHz, value in fraction
    225     ! has to be update once we have information from wave model (discussion with G. Madec)
    226      
    227     ! --------------------------------------------------------------------
    228     !  *- OSA estimation
    229     ! --------------------------------------------------------------------
    230     ! partitionning direct and diffuse albedo
    231     !
    232 
    233     ! excluding diffuse albedo ZRW on ZDIR_ALB
    234       ZDIR_ALB(JI,JWL) =  XFRWL(JWL) *((1.-ZFWC) * (ZR11+ZRW) +ZFWC*XRWC(JWL))
    235       ZSCA_ALB(JI,JWL) =  XFRWL(JWL) *((1.-ZFWC) * (ZRDF+ZRWDF) + ZFWC*XRWC(JWL))
    236 
    237     !  print*,ji,ZFWC,ZDIR_ALB(JI,JWL),ZSCA_ALB(JI,JWL),pwind(ji)
    238     ENDDO ! end of the loop over grid points
     202  CASE(4)
     203    IF (JWL.LE.49) THEN       ! from 200  to 680 nm
     204     inu=1
     205    ELSE IF (JWL.LE.99) THEN  ! from 690  to 1180 nm
     206     inu=2
     207    ELSE IF (JWL.LE.218) THEN ! from 1190 to 2370 nm
     208     inu=3
     209    ELSE                      ! from 2380 to 4000 nm
     210     inu=4
     211    ENDIF
     212  CASE(6)
     213    IF (JWL.LE.5) THEN        ! from 200  to 240 nm
     214     inu=1
     215    ELSE IF (JWL.LE.24) THEN  ! from 250  to 430 nm
     216     inu=2
     217    ELSE IF (JWL.LE.49) THEN  ! from 440  to 680 nm
     218     inu=3
     219    ELSE IF (JWL.LE.99) THEN  ! from 690  to 1180 nm
     220     inu=4
     221    ELSE IF (JWL.LE.218) THEN ! from 1190 to 2370 nm
     222     inu=5
     223    ELSE                      ! from 2380 to 4000 nm
     224     inu=6
     225    ENDIF
     226  END SELECT
     227
     228  ! partitionning direct and diffuse albedo
     229  ! excluding diffuse albedo ZRW on ZDIR_ALB
     230
     231  !--direct
     232  alb_dir_new(1:knon,inu)=alb_dir_new(1:knon,inu) + &
     233                          ( XFRWL(JWL) * ((1.-ZFWC(1:knon)) * (ZR11(1:knon)+ZRW(1:knon))   + ZFWC(1:knon)*XRWC(JWL)) )/SFRWL(inu)
     234  !--diffuse
     235  alb_dif_new(1:knon,inu)=alb_dif_new(1:knon,inu) + &
     236                          ( XFRWL(JWL) * ((1.-ZFWC(1:knon)) * (ZRDF(1:knon)+ZRWDF(1:knon)) + ZFWC(1:knon)*XRWC(JWL)) )/SFRWL(inu)
    239237
    240238ENDDO ! ending loop over wavelengths
    241239
    242 
    243 ! integral for each nsw band
    244 
    245 select case(nsw)
    246 case(2)
    247   do ji=1,knon
    248    alb_dir_new(ji,1)=sum(zdir_alb(ji,1:49))/SFRWL(1) ! from 200nm to 680nm
    249    alb_dir_new(ji,2)=sum(zdir_alb(ji,50:381))/SFRWL(2) ! from 690nm to 4000 nm
    250 
    251    alb_dif_new(ji,1)=sum(zsca_alb(ji,1:49))/SFRWL(1) ! from 200nm to 680nm
    252    alb_dif_new(ji,2)=sum(zsca_alb(ji,50:381))/SFRWL(2) ! from 690nm to 4000 nm
    253   enddo
    254 case(4)
    255   do ji=1,knon
    256    alb_dir_new(ji,1)=sum(zdir_alb(ji,1:49))/SFRWL(1) ! from 200nm to 680nm
    257    alb_dir_new(ji,2)=sum(zdir_alb(ji,50:99))/SFRWL(2) ! from 690nm to 1180 nm
    258    alb_dir_new(ji,3)=sum(zdir_alb(ji,100:218))/SFRWL(3) ! from 1190nm to 2370 nm
    259    alb_dir_new(ji,4)=sum(zdir_alb(ji,219:381))/SFRWL(4) ! from 2380nm to 4000 nm
    260 
    261    alb_dif_new(ji,1)=sum(zsca_alb(ji,1:49))/SFRWL(1) ! from 200nm to 680nm
    262    alb_dif_new(ji,2)=sum(zsca_alb(ji,50:99))/SFRWL(2) ! from 690nm to 1180 nm
    263    alb_dif_new(ji,3)=sum(zsca_alb(ji,100:218))/SFRWL(3) ! from 1190nm to 2370 nm
    264    alb_dif_new(ji,4)=sum(zsca_alb(ji,219:381))/SFRWL(4) ! from 2380nm to 4000 nm
    265   enddo
    266 case(6)
    267  do ji=1,knon
    268    alb_dir_new(ji,1)=sum(zdir_alb(ji,1:5))/SFRWL(1) ! from 200nm to 240nm
    269    alb_dir_new(ji,2)=sum(zdir_alb(ji,6:24))/SFRWL(2) ! from 250nm to 430 nm
    270    alb_dir_new(ji,3)=sum(zdir_alb(ji,25:49))/SFRWL(3) ! from 440nm to 680 nm
    271    alb_dir_new(ji,4)=sum(zdir_alb(ji,50:99))/SFRWL(4) ! from 690nm to 1180 nm
    272    alb_dir_new(ji,5)=sum(zdir_alb(ji,100:218))/SFRWL(5) ! from 1190nm to 2370 nm
    273    alb_dir_new(ji,6)=sum(zdir_alb(ji,219:381))/SFRWL(6) ! from 2380nm to 4000 nm
    274 
    275    alb_dif_new(ji,1)=sum(zsca_alb(ji,1:5))/SFRWL(1) ! from 200nm to 240nm
    276    alb_dif_new(ji,2)=sum(zsca_alb(ji,6:24))/SFRWL(2) ! from 250nm to 430 nm
    277    alb_dif_new(ji,3)=sum(zsca_alb(ji,25:49))/SFRWL(3) ! from 440nm to 680 nm
    278    alb_dif_new(ji,4)=sum(zsca_alb(ji,50:99))/SFRWL(4) ! from 690nm to 1180 nm
    279    alb_dif_new(ji,5)=sum(zsca_alb(ji,100:218))/SFRWL(5) ! from 1190nm to 2370 nm
    280    alb_dif_new(ji,6)=sum(zsca_alb(ji,219:381))/SFRWL(6) ! from 2380nm to 4000 nm
    281  enddo
    282 end select
    283 
    284 
    285 
    286 END subroutine ocean_albedo
    287 
     240END SUBROUTINE ocean_albedo
Note: See TracChangeset for help on using the changeset viewer.