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 8316 for branches/2017/dev_r8183_ICEMODEL/NEMOGCM/NEMO/OPA_SRC/SBC/sbccpl.F90 – NEMO

Ignore:
Timestamp:
2017-07-11T14:05:05+02:00 (7 years ago)
Author:
clem
Message:

STEP2 (3): remove obsolete features

File:
1 edited

Legend:

Unmodified
Added
Removed
  • branches/2017/dev_r8183_ICEMODEL/NEMOGCM/NEMO/OPA_SRC/SBC/sbccpl.F90

    r8313 r8316  
    5454 
    5555   PUBLIC   sbc_cpl_init      ! routine called by sbcmod.F90 
    56    PUBLIC   sbc_cpl_rcv       ! routine called by sbc_ice_lim(_2).F90 
     56   PUBLIC   sbc_cpl_rcv       ! routine called by sbc_ice_lim.F90 
    5757   PUBLIC   sbc_cpl_snd       ! routine called by step.F90 
    58    PUBLIC   sbc_cpl_ice_tau   ! routine called by sbc_ice_lim(_2).F90 
    59    PUBLIC   sbc_cpl_ice_flx   ! routine called by sbc_ice_lim(_2).F90 
     58   PUBLIC   sbc_cpl_ice_tau   ! routine called by sbc_ice_lim.F90 
     59   PUBLIC   sbc_cpl_ice_flx   ! routine called by sbc_ice_lim.F90 
    6060   PUBLIC   sbc_cpl_alloc     ! routine called in sbcice_cice.F90 
    6161 
     
    500500      ! 
    501501      ! non solar sensitivity mandatory for LIM ice model 
    502       IF( TRIM( sn_rcv_dqnsdt%cldes ) == 'none' .AND. k_ice /= 0 .AND. k_ice /= 4 .AND. nn_components /= jp_iam_sas ) & 
     502      IF( TRIM( sn_rcv_dqnsdt%cldes ) == 'none' .AND. k_ice /= 0 .AND. k_ice /= 3 .AND. nn_components /= jp_iam_sas ) & 
    503503         CALL ctl_stop( 'sbc_cpl_init: sn_rcv_dqnsdt%cldes must be coupled in namsbc_cpl namelist' ) 
    504504      ! non solar sensitivity mandatory for mixed oce-ice solar radiation coupling technique 
     
    15241524    
    15251525 
    1526    SUBROUTINE sbc_cpl_ice_flx( p_frld, palbi, psst, pist ) 
     1526   SUBROUTINE sbc_cpl_ice_flx( picefr, palbi, psst, pist ) 
    15271527      !!---------------------------------------------------------------------- 
    15281528      !!             ***  ROUTINE sbc_cpl_ice_flx  *** 
     
    15741574      !!                   sprecip           solid precipitation over the ocean   
    15751575      !!---------------------------------------------------------------------- 
    1576       REAL(wp), INTENT(in   ), DIMENSION(:,:)   ::   p_frld     ! lead fraction                [0 to 1] 
     1576      REAL(wp), INTENT(in), DIMENSION(:,:)             ::   picefr     ! ice fraction                [0 to 1] 
    15771577      ! optional arguments, used only in 'mixed oce-ice' case 
    1578       REAL(wp), INTENT(in   ), DIMENSION(:,:,:), OPTIONAL ::   palbi      ! all skies ice albedo  
    1579       REAL(wp), INTENT(in   ), DIMENSION(:,:  ), OPTIONAL ::   psst       ! sea surface temperature     [Celsius] 
    1580       REAL(wp), INTENT(in   ), DIMENSION(:,:,:), OPTIONAL ::   pist       ! ice surface temperature     [Kelvin] 
     1578      REAL(wp), INTENT(in), DIMENSION(:,:,:), OPTIONAL ::   palbi      ! all skies ice albedo  
     1579      REAL(wp), INTENT(in), DIMENSION(:,:  ), OPTIONAL ::   psst       ! sea surface temperature     [Celsius] 
     1580      REAL(wp), INTENT(in), DIMENSION(:,:,:), OPTIONAL ::   pist       ! ice surface temperature     [Kelvin] 
    15811581      ! 
    15821582      INTEGER ::   jl         ! dummy loop index 
    1583       REAL(wp), POINTER, DIMENSION(:,:  ) ::   zcptn, zcptrain, zcptsnw, zicefr, zmsk, zsnw 
     1583      REAL(wp), POINTER, DIMENSION(:,:  ) ::   zcptn, zcptrain, zcptsnw, ziceld, zmsk, zsnw 
    15841584      REAL(wp), POINTER, DIMENSION(:,:  ) ::   zemp_tot, zemp_ice, zemp_oce, ztprecip, zsprecip, zevap_oce, zevap_ice, zdevap_ice 
    15851585      REAL(wp), POINTER, DIMENSION(:,:  ) ::   zqns_tot, zqns_oce, zqsr_tot, zqsr_oce, zqprec_ice, zqemp_oce, zqemp_ice 
     
    15891589      IF( nn_timing == 1 )  CALL timing_start('sbc_cpl_ice_flx') 
    15901590      ! 
    1591       CALL wrk_alloc( jpi,jpj,     zcptn, zcptrain, zcptsnw, zicefr, zmsk, zsnw ) 
     1591      CALL wrk_alloc( jpi,jpj,     zcptn, zcptrain, zcptsnw, ziceld, zmsk, zsnw ) 
    15921592      CALL wrk_alloc( jpi,jpj,     zemp_tot, zemp_ice, zemp_oce, ztprecip, zsprecip, zevap_oce, zevap_ice, zdevap_ice ) 
    15931593      CALL wrk_alloc( jpi,jpj,     zqns_tot, zqns_oce, zqsr_tot, zqsr_oce, zqprec_ice, zqemp_oce, zqemp_ice ) 
     
    15951595 
    15961596      IF( ln_mixcpl )   zmsk(:,:) = 1. - xcplmask(:,:,0) 
    1597       zicefr(:,:) = 1.- p_frld(:,:) 
     1597      ziceld(:,:) = 1. - picefr(:,:) 
    15981598      zcptn(:,:) = rcp * sst_m(:,:) 
    15991599      ! 
     
    16111611         ztprecip(:,:) =   frcv(jpr_rain)%z3(:,:,1) + zsprecip(:,:)  ! May need to ensure positive here 
    16121612         zemp_tot(:,:) =   frcv(jpr_tevp)%z3(:,:,1) - ztprecip(:,:) 
    1613          zemp_ice(:,:) = ( frcv(jpr_ievp)%z3(:,:,1) - frcv(jpr_snow)%z3(:,:,1) ) * zicefr(:,:) 
     1613         zemp_ice(:,:) = ( frcv(jpr_ievp)%z3(:,:,1) - frcv(jpr_snow)%z3(:,:,1) ) * picefr(:,:) 
    16141614      CASE( 'oce and ice'   )   ! received fields: jpr_sbpr, jpr_semp, jpr_oemp, jpr_ievp 
    1615          zemp_tot(:,:) = p_frld(:,:) * frcv(jpr_oemp)%z3(:,:,1) + zicefr(:,:) * frcv(jpr_sbpr)%z3(:,:,1) 
    1616          zemp_ice(:,:) = frcv(jpr_semp)%z3(:,:,1) * zicefr(:,:) 
     1615         zemp_tot(:,:) = ziceld(:,:) * frcv(jpr_oemp)%z3(:,:,1) + picefr(:,:) * frcv(jpr_sbpr)%z3(:,:,1) 
     1616         zemp_ice(:,:) = frcv(jpr_semp)%z3(:,:,1) * picefr(:,:) 
    16171617         zsprecip(:,:) = frcv(jpr_ievp)%z3(:,:,1) - frcv(jpr_semp)%z3(:,:,1) 
    16181618         ztprecip(:,:) = frcv(jpr_semp)%z3(:,:,1) - frcv(jpr_sbpr)%z3(:,:,1) + zsprecip(:,:) 
     
    16201620 
    16211621#if defined key_lim3 
    1622       ! zsnw = snow fraction over ice after wind blowing (=zicefr if no blowing) 
    1623       zsnw(:,:) = 0._wp  ;  CALL lim_thd_snwblow( p_frld, zsnw ) 
     1622      ! zsnw = snow fraction over ice after wind blowing (=picefr if no blowing) 
     1623      zsnw(:,:) = 0._wp  ;  CALL lim_thd_snwblow( ziceld, zsnw ) 
    16241624       
    16251625      ! --- evaporation minus precipitation corrected (because of wind blowing on snow) --- ! 
    1626       zemp_ice(:,:) = zemp_ice(:,:) + zsprecip(:,:) * ( zicefr(:,:) - zsnw(:,:) )  ! emp_ice = A * sublimation - zsnw * sprecip 
     1626      zemp_ice(:,:) = zemp_ice(:,:) + zsprecip(:,:) * ( picefr(:,:) - zsnw(:,:) )  ! emp_ice = A * sublimation - zsnw * sprecip 
    16271627      zemp_oce(:,:) = zemp_tot(:,:) - zemp_ice(:,:)                                ! emp_oce = emp_tot - emp_ice 
    16281628 
    16291629      ! --- evaporation over ocean (used later for qemp) --- ! 
    1630       zevap_oce(:,:) = frcv(jpr_tevp)%z3(:,:,1) - frcv(jpr_ievp)%z3(:,:,1) * zicefr(:,:) 
     1630      zevap_oce(:,:) = frcv(jpr_tevp)%z3(:,:,1) - frcv(jpr_ievp)%z3(:,:,1) * picefr(:,:) 
    16311631 
    16321632      ! --- evaporation over ice (kg/m2/s) --- ! 
     
    16751675 
    16761676#else 
    1677       zsnw(:,:) = zicefr(:,:) 
     1677      zsnw(:,:) = picefr(:,:) 
    16781678      ! --- Continental fluxes --- ! 
    16791679      IF( srcv(jpr_rnf)%laction ) THEN   ! runoffs (included in emp later on) 
     
    17141714      IF( iom_use('snow_ao_cea') )  CALL iom_put( 'snow_ao_cea' , sprecip(:,:) * ( 1._wp - zsnw(:,:) )                  )  ! Snow over ice-free ocean  (cell average) 
    17151715      IF( iom_use('snow_ai_cea') )  CALL iom_put( 'snow_ai_cea' , sprecip(:,:) *           zsnw(:,:)                    )  ! Snow over sea-ice         (cell average) 
    1716       IF( iom_use('subl_ai_cea') )  CALL iom_put( 'subl_ai_cea' , frcv(jpr_ievp)%z3(:,:,1) * zicefr(:,:) * tmask(:,:,1) )  ! Sublimation over sea-ice (cell average) 
     1716      IF( iom_use('subl_ai_cea') )  CALL iom_put( 'subl_ai_cea' , frcv(jpr_ievp)%z3(:,:,1) * picefr(:,:) * tmask(:,:,1) )  ! Sublimation over sea-ice (cell average) 
    17171717      IF( iom_use('evap_ao_cea') )  CALL iom_put( 'evap_ao_cea' , ( frcv(jpr_tevp)%z3(:,:,1)  & 
    1718          &                                                        - frcv(jpr_ievp)%z3(:,:,1) * zicefr(:,:) ) * tmask(:,:,1) )  ! ice-free oce evap (cell average) 
     1718         &                                                        - frcv(jpr_ievp)%z3(:,:,1) * picefr(:,:) ) * tmask(:,:,1) )  ! ice-free oce evap (cell average) 
    17191719      ! note: runoff output is done in sbcrnf (which includes icebergs too) and iceshelf output is done in sbcisf 
    17201720      ! 
     
    17341734         ENDIF 
    17351735      CASE( 'oce and ice' )      ! the total flux is computed from ocean and ice fluxes 
    1736          zqns_tot(:,:) =  p_frld(:,:) * frcv(jpr_qnsoce)%z3(:,:,1) 
     1736         zqns_tot(:,:) =  ziceld(:,:) * frcv(jpr_qnsoce)%z3(:,:,1) 
    17371737         IF ( TRIM(sn_rcv_qns%clcat) == 'yes' ) THEN 
    17381738            DO jl=1,jpl 
     
    17411741            ENDDO 
    17421742         ELSE 
    1743             qns_tot(:,:) = qns_tot(:,:) + zicefr(:,:) * frcv(jpr_qnsice)%z3(:,:,1) 
     1743            qns_tot(:,:) = qns_tot(:,:) + picefr(:,:) * frcv(jpr_qnsice)%z3(:,:,1) 
    17441744            DO jl=1,jpl 
    1745                zqns_tot(:,:   ) = zqns_tot(:,:) + zicefr(:,:) * frcv(jpr_qnsice)%z3(:,:,1) 
     1745               zqns_tot(:,:   ) = zqns_tot(:,:) + picefr(:,:) * frcv(jpr_qnsice)%z3(:,:,1) 
    17461746               zqns_ice(:,:,jl) = frcv(jpr_qnsice)%z3(:,:,1) 
    17471747            ENDDO 
     
    17511751         zqns_tot(:,:  ) = frcv(jpr_qnsmix)%z3(:,:,1) 
    17521752         zqns_ice(:,:,1) = frcv(jpr_qnsmix)%z3(:,:,1)    & 
    1753             &            + frcv(jpr_dqnsdt)%z3(:,:,1) * ( pist(:,:,1) - ( (rt0 + psst(:,:  ) ) * p_frld(:,:)   & 
    1754             &                                           + pist(:,:,1) * zicefr(:,:) ) ) 
     1753            &            + frcv(jpr_dqnsdt)%z3(:,:,1) * ( pist(:,:,1) - ( (rt0 + psst(:,:  ) ) * ziceld(:,:)   & 
     1754            &                                           + pist(:,:,1) * picefr(:,:) ) ) 
    17551755      END SELECT 
    17561756      !                                      
     
    17631763#if defined key_lim3       
    17641764      ! --- non solar flux over ocean --- ! 
    1765       !         note: p_frld cannot be = 0 since we limit the ice concentration to amax 
     1765      !         note: ziceld cannot be = 0 since we limit the ice concentration to amax 
    17661766      zqns_oce = 0._wp 
    1767       WHERE( p_frld /= 0._wp )  zqns_oce(:,:) = ( zqns_tot(:,:) - SUM( a_i * zqns_ice, dim=3 ) ) / p_frld(:,:) 
     1767      WHERE( ziceld /= 0._wp )  zqns_oce(:,:) = ( zqns_tot(:,:) - SUM( a_i * zqns_ice, dim=3 ) ) / ziceld(:,:) 
    17681768 
    17691769      ! Heat content per unit mass of snow (J/kg) 
     
    17721772      ENDWHERE 
    17731773      ! Heat content per unit mass of rain (J/kg) 
    1774       zcptrain(:,:) = rcp * ( SUM( (tn_ice(:,:,:) - rt0) * a_i(:,:,:), dim=3 ) + sst_m(:,:) * p_frld(:,:) )  
     1774      zcptrain(:,:) = rcp * ( SUM( (tn_ice(:,:,:) - rt0) * a_i(:,:,:), dim=3 ) + sst_m(:,:) * ziceld(:,:) )  
    17751775 
    17761776      ! --- enthalpy of snow precip over ice in J/m3 (to be used in 1D-thermo) --- ! 
     
    17871787         &             +   zsprecip(:,:)                   * ( 1._wp - zsnw ) * ( zcptsnw (:,:) - lfus )   ! solid precip over ocean + snow melting 
    17881788      zqemp_ice(:,:) =     zsprecip(:,:)                   * zsnw             * ( zcptsnw (:,:) - lfus )   ! solid precip over ice (qevap_ice=0 since atm. does not take it into account) 
    1789 !!    zqemp_ice(:,:) = -   frcv(jpr_ievp)%z3(:,:,1)        * zicefr(:,:)      *   zcptsnw (:,:)   &        ! ice evap 
     1789!!    zqemp_ice(:,:) = -   frcv(jpr_ievp)%z3(:,:,1)        * picefr(:,:)      *   zcptsnw (:,:)   &        ! ice evap 
    17901790!!       &             +   zsprecip(:,:)                   * zsnw             * zqprec_ice(:,:) * r1_rhosn ! solid precip over ice 
    17911791       
     
    18201820      ! clem: this formulation is certainly wrong... but better than it was... 
    18211821      zqns_tot(:,:) = zqns_tot(:,:)                            &          ! zqns_tot update over free ocean with: 
    1822          &          - (  p_frld(:,:) * zsprecip(:,:) * lfus )  &          ! remove the latent heat flux of solid precip. melting 
     1822         &          - (  ziceld(:,:) * zsprecip(:,:) * lfus )  &          ! remove the latent heat flux of solid precip. melting 
    18231823         &          - (  zemp_tot(:,:)                         &          ! remove the heat content of mass flux (assumed to be at SST) 
    18241824         &             - zemp_ice(:,:) ) * zcptn(:,:)  
    18251825 
    18261826     IF( ln_mixcpl ) THEN 
    1827          qns_tot(:,:) = qns(:,:) * p_frld(:,:) + SUM( qns_ice(:,:,:) * a_i(:,:,:), dim=3 )   ! total flux from blk 
     1827         qns_tot(:,:) = qns(:,:) * ziceld(:,:) + SUM( qns_ice(:,:,:) * a_i(:,:,:), dim=3 )   ! total flux from blk 
    18281828         qns_tot(:,:) = qns_tot(:,:) * xcplmask(:,:,0) +  zqns_tot(:,:)* zmsk(:,:) 
    18291829         DO jl=1,jpl 
     
    18411841      IF( iom_use('hflx_snow_cea') ) CALL iom_put('hflx_snow_cea',  sprecip(:,:) * ( zcptsnw(:,:) - Lfus )                           ) ! heat flux from snow (cell average) 
    18421842      IF( iom_use('hflx_rain_cea') ) CALL iom_put('hflx_rain_cea',( tprecip(:,:) - sprecip(:,:) ) * zcptrain(:,:)                    ) ! heat flux from rain (cell average) 
    1843       IF( iom_use('hflx_evap_cea') ) CALL iom_put('hflx_evap_cea',(frcv(jpr_tevp)%z3(:,:,1) - frcv(jpr_ievp)%z3(:,:,1) * zicefr(:,:) & ! heat flux from from evap (cell average) 
     1843      IF( iom_use('hflx_evap_cea') ) CALL iom_put('hflx_evap_cea',(frcv(jpr_tevp)%z3(:,:,1) - frcv(jpr_ievp)%z3(:,:,1) * picefr(:,:) & ! heat flux from from evap (cell average) 
    18441844         &                                                        ) * zcptn(:,:) * tmask(:,:,1) ) 
    18451845      IF( iom_use('hflx_snow_ao_cea') ) CALL iom_put('hflx_snow_ao_cea',sprecip(:,:) * (zcptsnw(:,:) - Lfus) * (1._wp - zsnw(:,:))   ) ! heat flux from snow (over ocean) 
     
    18651865         zqsr_ice(:,:,1) = frcv(jpr_qsrice)%z3(:,:,1) 
    18661866      CASE( 'oce and ice' ) 
    1867          zqsr_tot(:,:  ) =  p_frld(:,:) * frcv(jpr_qsroce)%z3(:,:,1) 
     1867         zqsr_tot(:,:  ) =  ziceld(:,:) * frcv(jpr_qsroce)%z3(:,:,1) 
    18681868         IF ( TRIM(sn_rcv_qsr%clcat) == 'yes' ) THEN 
    18691869            DO jl=1,jpl 
     
    18721872            ENDDO 
    18731873         ELSE 
    1874             qsr_tot(:,:   ) = qsr_tot(:,:) + zicefr(:,:) * frcv(jpr_qsrice)%z3(:,:,1) 
     1874            qsr_tot(:,:   ) = qsr_tot(:,:) + picefr(:,:) * frcv(jpr_qsrice)%z3(:,:,1) 
    18751875            DO jl=1,jpl 
    1876                zqsr_tot(:,:   ) = zqsr_tot(:,:) + zicefr(:,:) * frcv(jpr_qsrice)%z3(:,:,1) 
     1876               zqsr_tot(:,:   ) = zqsr_tot(:,:) + picefr(:,:) * frcv(jpr_qsrice)%z3(:,:,1) 
    18771877               zqsr_ice(:,:,jl) = frcv(jpr_qsrice)%z3(:,:,1) 
    18781878            ENDDO 
     
    18841884!       ( see OASIS3 user guide, 5th edition, p39 ) 
    18851885         zqsr_ice(:,:,1) = frcv(jpr_qsrmix)%z3(:,:,1) * ( 1.- palbi(:,:,1) )   & 
    1886             &            / (  1.- ( albedo_oce_mix(:,:  ) * p_frld(:,:)       & 
    1887             &                     + palbi         (:,:,1) * zicefr(:,:) ) ) 
     1886            &            / (  1.- ( albedo_oce_mix(:,:  ) * ziceld(:,:)       & 
     1887            &                     + palbi         (:,:,1) * picefr(:,:) ) ) 
    18881888      END SELECT 
    18891889      IF( ln_dm2dc .AND. ln_cpl ) THEN   ! modify qsr to include the diurnal cycle 
     
    18961896#if defined key_lim3 
    18971897      ! --- solar flux over ocean --- ! 
    1898       !         note: p_frld cannot be = 0 since we limit the ice concentration to amax 
     1898      !         note: ziceld cannot be = 0 since we limit the ice concentration to amax 
    18991899      zqsr_oce = 0._wp 
    1900       WHERE( p_frld /= 0._wp )  zqsr_oce(:,:) = ( zqsr_tot(:,:) - SUM( a_i * zqsr_ice, dim=3 ) ) / p_frld(:,:) 
     1900      WHERE( ziceld /= 0._wp )  zqsr_oce(:,:) = ( zqsr_tot(:,:) - SUM( a_i * zqsr_ice, dim=3 ) ) / ziceld(:,:) 
    19011901 
    19021902      IF( ln_mixcpl ) THEN   ;   qsr_oce(:,:) = qsr_oce(:,:) * xcplmask(:,:,0) +  zqsr_oce(:,:)* zmsk(:,:) 
     
    19051905 
    19061906      IF( ln_mixcpl ) THEN 
    1907          qsr_tot(:,:) = qsr(:,:) * p_frld(:,:) + SUM( qsr_ice(:,:,:) * a_i(:,:,:), dim=3 )   ! total flux from blk 
     1907         qsr_tot(:,:) = qsr(:,:) * ziceld(:,:) + SUM( qsr_ice(:,:,:) * a_i(:,:,:), dim=3 )   ! total flux from blk 
    19081908         qsr_tot(:,:) = qsr_tot(:,:) * xcplmask(:,:,0) +  zqsr_tot(:,:)* zmsk(:,:) 
    19091909         DO jl=1,jpl 
     
    19521952      fr2_i0(:,:) = ( 0.82 * ( 1.0 - cldf_ice ) + 0.65 * cldf_ice ) 
    19531953 
    1954       CALL wrk_dealloc( jpi,jpj,     zcptn, zcptrain, zcptsnw, zicefr, zmsk, zsnw ) 
     1954      CALL wrk_dealloc( jpi,jpj,     zcptn, zcptrain, zcptsnw, ziceld, zmsk, zsnw ) 
    19551955      CALL wrk_dealloc( jpi,jpj,     zemp_tot, zemp_ice, zemp_oce, ztprecip, zsprecip, zevap_oce, zevap_ice, zdevap_ice ) 
    19561956      CALL wrk_dealloc( jpi,jpj,     zqns_tot, zqns_oce, zqsr_tot, zqsr_oce, zqprec_ice, zqemp_oce, zqemp_ice ) 
Note: See TracChangeset for help on using the changeset viewer.