New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
Changeset 14595 for NEMO/trunk/src/OCE – NEMO

Changeset 14595 for NEMO/trunk/src/OCE


Ignore:
Timestamp:
2021-03-05T23:36:50+01:00 (3 years ago)
Author:
clem
Message:

trunk: solve ticket #2627.

Location:
NEMO/trunk/src/OCE/SBC
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • NEMO/trunk/src/OCE/SBC/sbccpl.F90

    r14433 r14595  
    129129   INTEGER, PARAMETER ::   jpr_icb    = 61 
    130130   INTEGER, PARAMETER ::   jpr_ts_ice = 62   ! Sea ice surface temp 
     131   !!INTEGER, PARAMETER ::   jpr_qtrice = 63   ! Transmitted solar thru sea-ice 
    131132 
    132133   INTEGER, PARAMETER ::   jprcv      = 62   ! total number of fields received 
     
    202203      &             sn_rcv_wstrf, sn_rcv_wdrag, sn_rcv_charn, sn_rcv_taw, sn_rcv_bhd, sn_rcv_tusd, sn_rcv_tvsd 
    203204   !                                   ! Other namelist parameters 
     205!!   TYPE(FLD_C) ::   sn_rcv_qtrice 
    204206   INTEGER     ::   nn_cplmodel           ! Maximum number of models to/from which NEMO is potentialy sending/receiving data 
    205207   LOGICAL     ::   ln_usecplmask         !  use a coupling mask file to merge data received from several models 
     
    237239      !!             ***  FUNCTION sbc_cpl_alloc  *** 
    238240      !!---------------------------------------------------------------------- 
    239       INTEGER :: ierr(5) 
     241      INTEGER :: ierr(4) 
    240242      !!---------------------------------------------------------------------- 
    241243      ierr(:) = 0 
     
    247249#endif 
    248250      ALLOCATE( xcplmask(jpi,jpj,0:nn_cplmodel) , STAT=ierr(3) ) 
    249 #if defined key_si3 || defined key_cice 
    250       ALLOCATE( a_i_last_couple(jpi,jpj,jpl) , STAT=ierr(4) ) 
    251 #endif 
    252       ! 
    253       IF( .NOT. ln_apr_dyn ) ALLOCATE( ssh_ib(jpi,jpj), ssh_ibb(jpi,jpj), apr(jpi, jpj), STAT=ierr(5) ) 
     251      ! 
     252      IF( .NOT. ln_apr_dyn ) ALLOCATE( ssh_ib(jpi,jpj), ssh_ibb(jpi,jpj), apr(jpi, jpj), STAT=ierr(4) ) 
    254253 
    255254      sbc_cpl_alloc = MAXVAL( ierr ) 
     
    286285         &                  sn_rcv_charn , sn_rcv_taw   , sn_rcv_bhd  , sn_rcv_tusd  , sn_rcv_tvsd,    & 
    287286         &                  sn_rcv_wdrag , sn_rcv_qns   , sn_rcv_emp  , sn_rcv_rnf   , sn_rcv_cal  ,   & 
    288          &                  sn_rcv_iceflx, sn_rcv_co2   , sn_rcv_icb  , sn_rcv_isf   , sn_rcv_ts_ice 
     287         &                  sn_rcv_iceflx, sn_rcv_co2   , sn_rcv_icb  , sn_rcv_isf   , sn_rcv_ts_ice !!, sn_rcv_qtrice 
    289288 
    290289      !!--------------------------------------------------------------------- 
     
    327326         WRITE(numout,*)'      ice shelf                       = ', TRIM(sn_rcv_isf%cldes   ), ' (', TRIM(sn_rcv_isf%clcat   ), ')' 
    328327         WRITE(numout,*)'      sea ice heat fluxes             = ', TRIM(sn_rcv_iceflx%cldes), ' (', TRIM(sn_rcv_iceflx%clcat), ')' 
     328!!       WRITE(numout,*)'      transmitted solar thru sea-ice  = ', TRIM(sn_rcv_qtrice%cldes), ' (', TRIM(sn_rcv_qtrice%clcat), ')' 
    329329         WRITE(numout,*)'      atm co2                         = ', TRIM(sn_rcv_co2%cldes   ), ' (', TRIM(sn_rcv_co2%clcat   ), ')' 
    330330         WRITE(numout,*)'      Sea ice surface skin temperature= ', TRIM(sn_rcv_ts_ice%cldes), ' (', TRIM(sn_rcv_ts_ice%clcat), ')' 
     
    602602      srcv(jpr_mslp)%clname = 'O_MSLP'     ;   IF( TRIM(sn_rcv_mslp%cldes  ) == 'coupled' )    srcv(jpr_mslp)%laction = .TRUE. 
    603603      ! 
    604       !                                                      ! ------------------------- ! 
    605       !                                                      !  ice topmelt and botmelt  ! 
    606       !                                                      ! ------------------------- ! 
     604      !                                                      ! --------------------------------- ! 
     605      !                                                      !  ice topmelt and conduction flux  !    
     606      !                                                      ! --------------------------------- ! 
    607607      srcv(jpr_topm )%clname = 'OTopMlt' 
    608608      srcv(jpr_botm )%clname = 'OBotMlt' 
     
    615615         srcv(jpr_topm:jpr_botm)%laction = .TRUE. 
    616616      ENDIF 
     617!!      !                                                      ! --------------------------- ! 
     618!!      !                                                      ! transmitted solar thru ice  !    
     619!!      !                                                      ! --------------------------- ! 
     620!!      srcv(jpr_qtrice)%clname = 'OQtr' 
     621!!      IF( TRIM(sn_rcv_qtrice%cldes) == 'coupled' ) THEN 
     622!!         IF ( TRIM( sn_rcv_qtrice%clcat ) == 'yes' ) THEN 
     623!!            srcv(jpr_qtrice)%nct = nn_cats_cpl 
     624!!         ELSE 
     625!!           CALL ctl_stop( 'sbc_cpl_init: sn_rcv_qtrice%clcat should always be set to yes currently' ) 
     626!!         ENDIF 
     627!!         srcv(jpr_qtrice)%laction = .TRUE. 
     628!!      ENDIF 
    617629      !                                                      ! ------------------------- ! 
    618630      !                                                      !    ice skin temperature   ! 
     
    888900      END SELECT 
    889901 
    890       ! Initialise ice fractions from last coupling time to zero (needed by Met-Office) 
    891 #if defined key_si3 || defined key_cice 
    892        a_i_last_couple(:,:,:) = 0._wp 
    893 #endif 
    894902      !                                                      ! ------------------------- ! 
    895903      !                                                      !      Ice Meltponds        ! 
     
    15891597      !! ** Action  :   return ptau_i, ptau_j, the stress over the ice 
    15901598      !!---------------------------------------------------------------------- 
    1591       REAL(wp), INTENT(out), DIMENSION(:,:) ::   p_taui   ! i- & j-components of atmos-ice stress [N/m2] 
    1592       REAL(wp), INTENT(out), DIMENSION(:,:) ::   p_tauj   ! at I-point (B-grid) or U & V-point (C-grid) 
     1599      REAL(wp), INTENT(inout), DIMENSION(:,:) ::   p_taui   ! i- & j-components of atmos-ice stress [N/m2] 
     1600      REAL(wp), INTENT(inout), DIMENSION(:,:) ::   p_tauj   ! at I-point (B-grid) or U & V-point (C-grid) 
    15931601      !! 
    15941602      INTEGER ::   ji, jj   ! dummy loop indices 
     
    15971605      REAL(wp), DIMENSION(jpi,jpj) ::   ztx, zty 
    15981606      !!---------------------------------------------------------------------- 
     1607      ! 
     1608#if defined key_si3 || defined key_cice 
    15991609      ! 
    16001610      IF( srcv(jpr_itx1)%laction ) THEN   ;   itx =  jpr_itx1 
     
    16711681      ENDIF 
    16721682      ! 
     1683#endif 
     1684      ! 
    16731685   END SUBROUTINE sbc_cpl_ice_tau 
    16741686 
    16751687 
    1676    SUBROUTINE sbc_cpl_ice_flx( picefr, palbi, psst, pist, phs, phi ) 
     1688   SUBROUTINE sbc_cpl_ice_flx( kt, picefr, palbi, psst, pist, phs, phi ) 
    16771689      !!---------------------------------------------------------------------- 
    16781690      !!             ***  ROUTINE sbc_cpl_ice_flx  *** 
     
    17161728      !!                                                                      are provided but not included in emp here. Only runoff will 
    17171729      !!                                                                      be included in emp in other parts of NEMO code 
     1730      !! 
     1731      !! ** Note : In case of the ice-atm coupling with conduction fluxes (such as Jules interface for the Met-Office), 
     1732      !!              qsr_ice and qns_ice are not provided and they are not supposed to be used in the ice code. 
     1733      !!              However, by precaution we also "fake" qns_ice and qsr_ice this way: 
     1734      !!              qns_ice = qml_ice + qcn_ice ?? 
     1735      !!              qsr_ice = qtr_ice_top ?? 
     1736      !! 
    17181737      !! ** Action  :   update at each nf_ice time step: 
    17191738      !!                   qns_tot, qsr_tot  non-solar and solar total heat fluxes 
     
    17241743      !!                   sprecip           solid precipitation over the ocean 
    17251744      !!---------------------------------------------------------------------- 
     1745      INTEGER,  INTENT(in)                                ::   kt         ! ocean model time step index (only for a_i_last_couple) 
    17261746      REAL(wp), INTENT(in)   , DIMENSION(:,:)             ::   picefr     ! ice fraction                [0 to 1] 
    17271747      !                                                   !!           ! optional arguments, used only in 'mixed oce-ice' case or for Met-Office coupling 
     
    17401760      REAL(wp), DIMENSION(jpi,jpj)     ::   ztri 
    17411761      !!---------------------------------------------------------------------- 
     1762      ! 
     1763#if defined key_si3 || defined key_cice 
     1764      ! 
     1765      IF( kt == nit000 ) THEN 
     1766         ! allocate ice fractions from last coupling time here and not in sbc_cpl_init because of jpl 
     1767         IF( .NOT.ALLOCATED(a_i_last_couple) )   ALLOCATE( a_i_last_couple(jpi,jpj,jpl) ) 
     1768         ! initialize to a_i for the 1st time step 
     1769         a_i_last_couple(:,:,:) = a_i(:,:,:) 
     1770      ENDIF 
    17421771      ! 
    17431772      IF( ln_mixcpl )   zmsk(:,:) = 1. - xcplmask(:,:,0) 
     
    17671796         CALL ctl_stop('STOP', 'sbccpl/sbc_cpl_ice_flx: some fields are not defined. Change sn_rcv_emp value in namelist namsbc_cpl') 
    17681797      END SELECT 
    1769  
    1770 #if defined key_si3 
    17711798 
    17721799      ! --- evaporation over ice (kg/m2/s) --- ! 
     
    18601887      ENDIF 
    18611888 
    1862 #else 
    1863       zsnw(:,:) = picefr(:,:) 
    1864       ! --- Continental fluxes --- ! 
    1865       IF( srcv(jpr_rnf)%laction ) THEN   ! runoffs (included in emp later on) 
    1866          rnf(:,:) = frcv(jpr_rnf)%z3(:,:,1) 
    1867       ENDIF 
    1868       IF( srcv(jpr_cal)%laction ) THEN   ! calving (put in emp_tot) 
    1869          zemp_tot(:,:) = zemp_tot(:,:) - frcv(jpr_cal)%z3(:,:,1) 
    1870       ENDIF 
    1871       IF( srcv(jpr_icb)%laction ) THEN   ! iceberg added to runoffs 
    1872          fwficb(:,:) = frcv(jpr_icb)%z3(:,:,1) 
    1873          rnf(:,:)    = rnf(:,:) + fwficb(:,:) 
    1874       ENDIF 
    1875       IF( srcv(jpr_isf)%laction ) THEN   ! iceshelf (fwfisf <0 mean melting) 
    1876         fwfisf_oasis(:,:) = - frcv(jpr_isf)%z3(:,:,1) 
    1877       ENDIF 
    1878       ! 
    1879       IF( ln_mixcpl ) THEN 
    1880          emp_tot(:,:) = emp_tot(:,:) * xcplmask(:,:,0) + zemp_tot(:,:) * zmsk(:,:) 
    1881          emp_ice(:,:) = emp_ice(:,:) * xcplmask(:,:,0) + zemp_ice(:,:) * zmsk(:,:) 
    1882          sprecip(:,:) = sprecip(:,:) * xcplmask(:,:,0) + zsprecip(:,:) * zmsk(:,:) 
    1883          tprecip(:,:) = tprecip(:,:) * xcplmask(:,:,0) + ztprecip(:,:) * zmsk(:,:) 
    1884       ELSE 
    1885          emp_tot(:,:) =                                  zemp_tot(:,:) 
    1886          emp_ice(:,:) =                                  zemp_ice(:,:) 
    1887          sprecip(:,:) =                                  zsprecip(:,:) 
    1888          tprecip(:,:) =                                  ztprecip(:,:) 
    1889       ENDIF 
    1890       ! 
    1891 #endif 
    1892  
     1889!! for CICE ?? 
     1890!!$      zsnw(:,:) = picefr(:,:) 
     1891!!$      ! --- Continental fluxes --- ! 
     1892!!$      IF( srcv(jpr_rnf)%laction ) THEN   ! runoffs (included in emp later on) 
     1893!!$         rnf(:,:) = frcv(jpr_rnf)%z3(:,:,1) 
     1894!!$      ENDIF 
     1895!!$      IF( srcv(jpr_cal)%laction ) THEN   ! calving (put in emp_tot) 
     1896!!$         zemp_tot(:,:) = zemp_tot(:,:) - frcv(jpr_cal)%z3(:,:,1) 
     1897!!$      ENDIF 
     1898!!$      IF( srcv(jpr_icb)%laction ) THEN   ! iceberg added to runoffs 
     1899!!$         fwficb(:,:) = frcv(jpr_icb)%z3(:,:,1) 
     1900!!$         rnf(:,:)    = rnf(:,:) + fwficb(:,:) 
     1901!!$      ENDIF 
     1902!!$      IF( srcv(jpr_isf)%laction ) THEN   ! iceshelf (fwfisf <0 mean melting) 
     1903!!$        fwfisf_oasis(:,:) = - frcv(jpr_isf)%z3(:,:,1) 
     1904!!$      ENDIF 
     1905!!$      ! 
     1906!!$      IF( ln_mixcpl ) THEN 
     1907!!$         emp_tot(:,:) = emp_tot(:,:) * xcplmask(:,:,0) + zemp_tot(:,:) * zmsk(:,:) 
     1908!!$         emp_ice(:,:) = emp_ice(:,:) * xcplmask(:,:,0) + zemp_ice(:,:) * zmsk(:,:) 
     1909!!$         sprecip(:,:) = sprecip(:,:) * xcplmask(:,:,0) + zsprecip(:,:) * zmsk(:,:) 
     1910!!$         tprecip(:,:) = tprecip(:,:) * xcplmask(:,:,0) + ztprecip(:,:) * zmsk(:,:) 
     1911!!$      ELSE 
     1912!!$         emp_tot(:,:) =                                  zemp_tot(:,:) 
     1913!!$         emp_ice(:,:) =                                  zemp_ice(:,:) 
     1914!!$         sprecip(:,:) =                                  zsprecip(:,:) 
     1915!!$         tprecip(:,:) =                                  ztprecip(:,:) 
     1916!!$      ENDIF 
     1917      ! 
    18931918      ! outputs 
    1894 !!      IF( srcv(jpr_rnf)%laction )   CALL iom_put( 'runoffs' , rnf(:,:) * tmask(:,:,1)                                 )  ! runoff 
    1895 !!      IF( srcv(jpr_isf)%laction )   CALL iom_put( 'iceshelf_cea', -fwfisf(:,:) * tmask(:,:,1)                         )  ! iceshelf 
    18961919      IF( srcv(jpr_cal)%laction )    CALL iom_put( 'calving_cea' , frcv(jpr_cal)%z3(:,:,1) * tmask(:,:,1)                )  ! calving 
    18971920      IF( srcv(jpr_icb)%laction )    CALL iom_put( 'iceberg_cea' , frcv(jpr_icb)%z3(:,:,1) * tmask(:,:,1)                )  ! icebergs 
     
    19061929         &                                                         - frcv(jpr_ievp)%z3(:,:,1) * picefr(:,:) ) * tmask(:,:,1) ) ! ice-free oce evap (cell average) 
    19071930      ! note: runoff output is done in sbcrnf (which includes icebergs too) and iceshelf output is done in sbcisf 
     1931!!      IF( srcv(jpr_rnf)%laction )   CALL iom_put( 'runoffs' , rnf(:,:) * tmask(:,:,1)                                 )  ! runoff 
     1932!!      IF( srcv(jpr_isf)%laction )   CALL iom_put( 'iceshelf_cea', -fwfisf(:,:) * tmask(:,:,1)                         )  ! iceshelf 
     1933      ! 
     1934      !                                                      ! ========================= ! 
     1935      SELECT CASE( TRIM( sn_rcv_iceflx%cldes ) )             !  ice topmelt and botmelt  ! 
     1936      !                                                      ! ========================= ! 
     1937      CASE ('coupled') 
     1938         IF (ln_scale_ice_flux) THEN 
     1939            WHERE( a_i(:,:,:) > 1.e-10_wp ) 
     1940               qml_ice(:,:,:) = frcv(jpr_topm)%z3(:,:,:) * a_i_last_couple(:,:,:) / a_i(:,:,:) 
     1941               qcn_ice(:,:,:) = frcv(jpr_botm)%z3(:,:,:) * a_i_last_couple(:,:,:) / a_i(:,:,:) 
     1942            ELSEWHERE 
     1943               qml_ice(:,:,:) = 0.0_wp 
     1944               qcn_ice(:,:,:) = 0.0_wp 
     1945            END WHERE 
     1946         ELSE 
     1947            qml_ice(:,:,:) = frcv(jpr_topm)%z3(:,:,:) 
     1948            qcn_ice(:,:,:) = frcv(jpr_botm)%z3(:,:,:) 
     1949         ENDIF 
     1950      END SELECT 
    19081951      ! 
    19091952      !                                                      ! ========================= ! 
     
    19111954      !                                                      ! ========================= ! 
    19121955      CASE( 'oce only' )         ! the required field is directly provided 
    1913          zqns_tot(:,:) = frcv(jpr_qnsoce)%z3(:,:,1) 
    1914          ! For Met Office sea ice non-solar fluxes are already delt with by JULES so setting to zero 
    1915          ! here so the only flux is the ocean only one. 
    1916          zqns_ice(:,:,:) = 0._wp 
     1956         ! Get the sea ice non solar heat flux from conductive, melting and sublimation fluxes 
     1957         IF( TRIM(sn_rcv_iceflx%cldes) == 'coupled' ) THEN 
     1958            zqns_ice(:,:,:) = qml_ice(:,:,:) + qcn_ice(:,:,:) 
     1959         ELSE 
     1960            zqns_ice(:,:,:) = 0._wp 
     1961         ENDIF 
     1962         ! Calculate the total non solar heat flux. The ocean only non solar heat flux (zqns_oce) will be recalculated after this CASE 
     1963         ! statement to be consistent with other coupling methods even though .zqns_oce = frcv(jpr_qnsoce)%z3(:,:,1) 
     1964         zqns_tot(:,:) = frcv(jpr_qnsoce)%z3(:,:,1) + SUM( zqns_ice(:,:,:) * a_i(:,:,:), dim=3 ) 
    19171965      CASE( 'conservative' )     ! the required fields are directly provided 
    19181966         zqns_tot(:,:) = frcv(jpr_qnsmix)%z3(:,:,1) 
     
    19612009      IF( srcv(jpr_icb)%laction )   zqns_tot(:,:) = zqns_tot(:,:) - frcv(jpr_icb)%z3(:,:,1) * rLfus  ! remove latent heat of iceberg melting 
    19622010 
    1963 #if defined key_si3 
    19642011      ! --- non solar flux over ocean --- ! 
    19652012      !         note: ziceld cannot be = 0 since we limit the ice concentration to amax 
     
    20142061      ENDIF 
    20152062 
    2016 #else 
    2017       zcptsnw (:,:) = zcptn(:,:) 
    2018       zcptrain(:,:) = zcptn(:,:) 
    2019  
    2020       ! clem: this formulation is certainly wrong... but better than it was... 
    2021       zqns_tot(:,:) = zqns_tot(:,:)                             &          ! zqns_tot update over free ocean with: 
    2022          &          - (  ziceld(:,:) * zsprecip(:,:) * rLfus )  &          ! remove the latent heat flux of solid precip. melting 
    2023          &          - (  zemp_tot(:,:)                          &          ! remove the heat content of mass flux (assumed to be at SST) 
    2024          &             - zemp_ice(:,:) ) * zcptn(:,:) 
    2025  
    2026      IF( ln_mixcpl ) THEN 
    2027          qns_tot(:,:) = qns(:,:) * ziceld(:,:) + SUM( qns_ice(:,:,:) * a_i(:,:,:), dim=3 )   ! total flux from blk 
    2028          qns_tot(:,:) = qns_tot(:,:) * xcplmask(:,:,0) +  zqns_tot(:,:)* zmsk(:,:) 
    2029          DO jl=1,jpl 
    2030             qns_ice(:,:,jl) = qns_ice(:,:,jl) * xcplmask(:,:,0) +  zqns_ice(:,:,jl)* zmsk(:,:) 
    2031          ENDDO 
    2032       ELSE 
    2033          qns_tot(:,:  ) = zqns_tot(:,:  ) 
    2034          qns_ice(:,:,:) = zqns_ice(:,:,:) 
    2035       ENDIF 
    2036  
    2037 #endif 
     2063!! for CICE ?? 
     2064!!$      ! --- non solar flux over ocean --- ! 
     2065!!$      zcptsnw (:,:) = zcptn(:,:) 
     2066!!$      zcptrain(:,:) = zcptn(:,:) 
     2067!!$ 
     2068!!$      ! clem: this formulation is certainly wrong... but better than it was... 
     2069!!$      zqns_tot(:,:) = zqns_tot(:,:)                             &          ! zqns_tot update over free ocean with: 
     2070!!$         &          - (  ziceld(:,:) * zsprecip(:,:) * rLfus )  &          ! remove the latent heat flux of solid precip. melting 
     2071!!$         &          - (  zemp_tot(:,:)                          &          ! remove the heat content of mass flux (assumed to be at SST) 
     2072!!$         &             - zemp_ice(:,:) ) * zcptn(:,:) 
     2073!!$ 
     2074!!$     IF( ln_mixcpl ) THEN 
     2075!!$         qns_tot(:,:) = qns(:,:) * ziceld(:,:) + SUM( qns_ice(:,:,:) * a_i(:,:,:), dim=3 )   ! total flux from blk 
     2076!!$         qns_tot(:,:) = qns_tot(:,:) * xcplmask(:,:,0) +  zqns_tot(:,:)* zmsk(:,:) 
     2077!!$         DO jl=1,jpl 
     2078!!$            qns_ice(:,:,jl) = qns_ice(:,:,jl) * xcplmask(:,:,0) +  zqns_ice(:,:,jl)* zmsk(:,:) 
     2079!!$         ENDDO 
     2080!!$      ELSE 
     2081!!$         qns_tot(:,:  ) = zqns_tot(:,:  ) 
     2082!!$         qns_ice(:,:,:) = zqns_ice(:,:,:) 
     2083!!$      ENDIF 
     2084 
    20382085      ! outputs 
    20392086      IF ( srcv(jpr_cal)%laction ) CALL iom_put('hflx_cal_cea' , - frcv(jpr_cal)%z3(:,:,1) * rLfus ) ! latent heat from calving 
     
    20562103      ! 
    20572104      !                                                      ! ========================= ! 
     2105      SELECT CASE( TRIM( sn_rcv_dqnsdt%cldes ) )             !          d(qns)/dt        ! 
     2106      !                                                      ! ========================= ! 
     2107      CASE ('coupled') 
     2108         IF( TRIM(sn_rcv_dqnsdt%clcat) == 'yes' ) THEN 
     2109            zdqns_ice(:,:,1:jpl) = frcv(jpr_dqnsdt)%z3(:,:,1:jpl) 
     2110         ELSE 
     2111            ! Set all category values equal for the moment 
     2112            DO jl=1,jpl 
     2113               zdqns_ice(:,:,jl) = frcv(jpr_dqnsdt)%z3(:,:,1) 
     2114            ENDDO 
     2115         ENDIF 
     2116      CASE( 'none' ) 
     2117         zdqns_ice(:,:,:) = 0._wp 
     2118      END SELECT 
     2119 
     2120      IF( ln_mixcpl ) THEN 
     2121         DO jl=1,jpl 
     2122            dqns_ice(:,:,jl) = dqns_ice(:,:,jl) * xcplmask(:,:,0) + zdqns_ice(:,:,jl) * zmsk(:,:) 
     2123         ENDDO 
     2124      ELSE 
     2125         dqns_ice(:,:,:) = zdqns_ice(:,:,:) 
     2126      ENDIF 
     2127      !                                                      ! ========================= ! 
    20582128      SELECT CASE( TRIM( sn_rcv_qsr%cldes ) )                !      solar heat fluxes    !   (qsr) 
    20592129      !                                                      ! ========================= ! 
    20602130      CASE( 'oce only' ) 
    20612131         zqsr_tot(:,:  ) = MAX( 0._wp , frcv(jpr_qsroce)%z3(:,:,1) ) 
    2062          ! For Met Office sea ice solar fluxes are already delt with by JULES so setting to zero 
    2063          ! here so the only flux is the ocean only one. 
     2132         ! For the Met Office the only sea ice solar flux is the transmitted qsr which is added onto zqsr_ice 
     2133         ! further down. Therefore start zqsr_ice off at zero. 
    20642134         zqsr_ice(:,:,:) = 0._wp 
    20652135      CASE( 'conservative' ) 
     
    21142184         END DO 
    21152185      ENDIF 
    2116  
    2117 #if defined key_si3 
    2118       ! --- solar flux over ocean --- ! 
    2119       !         note: ziceld cannot be = 0 since we limit the ice concentration to amax 
    2120       zqsr_oce = 0._wp 
    2121       WHERE( ziceld /= 0._wp )  zqsr_oce(:,:) = ( zqsr_tot(:,:) - SUM( a_i * zqsr_ice, dim=3 ) ) / ziceld(:,:) 
    2122  
    2123       IF( ln_mixcpl ) THEN   ;   qsr_oce(:,:) = qsr_oce(:,:) * xcplmask(:,:,0) +  zqsr_oce(:,:)* zmsk(:,:) 
    2124       ELSE                   ;   qsr_oce(:,:) = zqsr_oce(:,:)   ;   ENDIF 
    2125 #endif 
    2126  
    2127       IF( ln_mixcpl ) THEN 
    2128          qsr_tot(:,:) = qsr(:,:) * ziceld(:,:) + SUM( qsr_ice(:,:,:) * a_i(:,:,:), dim=3 )   ! total flux from blk 
    2129          qsr_tot(:,:) = qsr_tot(:,:) * xcplmask(:,:,0) +  zqsr_tot(:,:)* zmsk(:,:) 
    2130          DO jl = 1, jpl 
    2131             qsr_ice(:,:,jl) = qsr_ice(:,:,jl) * xcplmask(:,:,0) +  zqsr_ice(:,:,jl)* zmsk(:,:) 
    2132          END DO 
    2133       ELSE 
    2134          qsr_tot(:,:  ) = zqsr_tot(:,:  ) 
    2135          qsr_ice(:,:,:) = zqsr_ice(:,:,:) 
    2136       ENDIF 
    2137  
    2138       !                                                      ! ========================= ! 
    2139       SELECT CASE( TRIM( sn_rcv_dqnsdt%cldes ) )             !          d(qns)/dt        ! 
    2140       !                                                      ! ========================= ! 
    2141       CASE ('coupled') 
    2142          IF( TRIM(sn_rcv_dqnsdt%clcat) == 'yes' ) THEN 
    2143             zdqns_ice(:,:,1:jpl) = frcv(jpr_dqnsdt)%z3(:,:,1:jpl) 
    2144          ELSE 
    2145             ! Set all category values equal for the moment 
    2146             DO jl=1,jpl 
    2147                zdqns_ice(:,:,jl) = frcv(jpr_dqnsdt)%z3(:,:,1) 
    2148             ENDDO 
    2149          ENDIF 
    2150       CASE( 'none' ) 
    2151          zdqns_ice(:,:,:) = 0._wp 
    2152       END SELECT 
    2153  
    2154       IF( ln_mixcpl ) THEN 
    2155          DO jl=1,jpl 
    2156             dqns_ice(:,:,jl) = dqns_ice(:,:,jl) * xcplmask(:,:,0) + zdqns_ice(:,:,jl) * zmsk(:,:) 
    2157          ENDDO 
    2158       ELSE 
    2159          dqns_ice(:,:,:) = zdqns_ice(:,:,:) 
    2160       ENDIF 
    2161  
    2162 #if defined key_si3 
    2163       !                                                      ! ========================= ! 
    2164       SELECT CASE( TRIM( sn_rcv_iceflx%cldes ) )             !  ice topmelt and botmelt  ! 
    2165       !                                                      ! ========================= ! 
    2166       CASE ('coupled') 
    2167          IF (ln_scale_ice_flux) THEN 
    2168             WHERE( a_i(:,:,:) > 1.e-10_wp ) 
    2169                qml_ice(:,:,:) = frcv(jpr_topm)%z3(:,:,:) * a_i_last_couple(:,:,:) / a_i(:,:,:) 
    2170                qcn_ice(:,:,:) = frcv(jpr_botm)%z3(:,:,:) * a_i_last_couple(:,:,:) / a_i(:,:,:) 
    2171             ELSEWHERE 
    2172                qml_ice(:,:,:) = 0.0_wp 
    2173                qcn_ice(:,:,:) = 0.0_wp 
    2174             END WHERE 
    2175          ELSE 
    2176             qml_ice(:,:,:) = frcv(jpr_topm)%z3(:,:,:) 
    2177             qcn_ice(:,:,:) = frcv(jpr_botm)%z3(:,:,:) 
    2178          ENDIF 
    2179       END SELECT 
    21802186      !                                                      ! ========================= ! 
    21812187      !                                                      !      Transmitted Qsr      !   [W/m2] 
     
    22092215      ELSEIF( ln_cndflx .AND. .NOT.ln_cndemulate ) THEN      !==  conduction flux as surface forcing  ==! 
    22102216         ! 
    2211          !          ! ===> here we must receive the qtr_ice_top array from the coupler 
    2212          !                 for now just assume zero (fully opaque ice) 
     2217!!         SELECT CASE( TRIM( sn_rcv_qtrice%cldes ) ) 
     2218!!            ! 
     2219!!            !      ! ===> here we receive the qtr_ice_top array from the coupler 
     2220!!         CASE ('coupled') 
     2221!!            IF (ln_scale_ice_flux) THEN 
     2222!!               WHERE( a_i(:,:,:) > 1.e-10_wp ) 
     2223!!                  zqtr_ice_top(:,:,:) = frcv(jpr_qtrice)%z3(:,:,:) * a_i_last_couple(:,:,:) / a_i(:,:,:) 
     2224!!               ELSEWHERE 
     2225!!                  zqtr_ice_top(:,:,:) = 0.0_wp 
     2226!!               ENDWHERE 
     2227!!            ELSE 
     2228!!               zqtr_ice_top(:,:,:) = frcv(jpr_qtrice)%z3(:,:,:) 
     2229!!            ENDIF 
     2230!!            
     2231!!            ! Add retrieved transmitted solar radiation onto the ice and total solar radiation 
     2232!!            zqsr_ice(:,:,:) = zqsr_ice(:,:,:) + zqtr_ice_top(:,:,:) 
     2233!!            zqsr_tot(:,:)   = zqsr_tot(:,:) + SUM( zqtr_ice_top(:,:,:) * a_i(:,:,:), dim=3 ) 
     2234!!             
     2235!!            !      if we are not getting this data from the coupler then assume zero (fully opaque ice) 
     2236!!         CASE ('none') 
    22132237         zqtr_ice_top(:,:,:) = 0._wp 
    2214          ! 
    2215       ENDIF 
    2216       ! 
     2238!!         END SELECT 
     2239            ! 
     2240      ENDIF 
     2241 
    22172242      IF( ln_mixcpl ) THEN 
    2218          DO jl=1,jpl 
     2243         qsr_tot(:,:) = qsr(:,:) * ziceld(:,:) + SUM( qsr_ice(:,:,:) * a_i(:,:,:), dim=3 )   ! total flux from blk 
     2244         qsr_tot(:,:) = qsr_tot(:,:) * xcplmask(:,:,0) + zqsr_tot(:,:) * zmsk(:,:) 
     2245         DO jl = 1, jpl 
     2246            qsr_ice    (:,:,jl) = qsr_ice    (:,:,jl) * xcplmask(:,:,0) + zqsr_ice    (:,:,jl) * zmsk(:,:) 
    22192247            qtr_ice_top(:,:,jl) = qtr_ice_top(:,:,jl) * xcplmask(:,:,0) + zqtr_ice_top(:,:,jl) * zmsk(:,:) 
    2220          ENDDO 
     2248         END DO 
    22212249      ELSE 
     2250         qsr_tot    (:,:  ) = zqsr_tot    (:,:  ) 
     2251         qsr_ice    (:,:,:) = zqsr_ice    (:,:,:) 
    22222252         qtr_ice_top(:,:,:) = zqtr_ice_top(:,:,:) 
    22232253      ENDIF 
     2254       
     2255      ! --- solar flux over ocean --- ! 
     2256      !         note: ziceld cannot be = 0 since we limit the ice concentration to amax 
     2257      zqsr_oce = 0._wp 
     2258      WHERE( ziceld /= 0._wp )  zqsr_oce(:,:) = ( zqsr_tot(:,:) - SUM( a_i * zqsr_ice, dim=3 ) ) / ziceld(:,:) 
     2259 
     2260      IF( ln_mixcpl ) THEN   ;   qsr_oce(:,:) = qsr_oce(:,:) * xcplmask(:,:,0) +  zqsr_oce(:,:)* zmsk(:,:) 
     2261      ELSE                   ;   qsr_oce(:,:) = zqsr_oce(:,:)   ;   ENDIF 
     2262 
    22242263      !                                                      ! ================== ! 
    22252264      !                                                      !   ice skin temp.   ! 
  • NEMO/trunk/src/OCE/SBC/sbcice_cice.F90

    r14433 r14595  
    139139            CALL cice_sbc_force(kt) 
    140140         ELSE IF( ksbc == jp_purecpl ) THEN 
    141             CALL sbc_cpl_ice_flx( fr_i ) 
     141            CALL sbc_cpl_ice_flx( kt, fr_i ) 
    142142         ENDIF 
    143143 
Note: See TracChangeset for help on using the changeset viewer.