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 13237 for NEMO/trunk/src/OCE/ISF – NEMO

Ignore:
Timestamp:
2020-07-03T11:12:53+02:00 (4 years ago)
Author:
smasson
Message:

trunk: Mid-year merge, merge back KERNEL-06_techene_e3

Location:
NEMO/trunk/src/OCE/ISF
Files:
8 edited

Legend:

Unmodified
Added
Removed
  • NEMO/trunk/src/OCE/ISF/isfcavgam.F90

    r12077 r13237  
    3030   PUBLIC   isfcav_gammats 
    3131 
     32#  include "domzgr_substitute.h90" 
    3233   !!---------------------------------------------------------------------- 
    3334   !! NEMO/OCE 4.0 , NEMO Consortium (2018) 
  • NEMO/trunk/src/OCE/ISF/isfcpl.F90

    r13226 r13237  
    1515   USE isfutils, ONLY : debug 
    1616   USE lib_mpp , ONLY: mpp_sum, mpp_max ! mpp routine 
     17#if ! defined key_qco 
    1718   USE domvvl  , ONLY: dom_vvl_zgr      ! vertical scale factor interpolation 
     19#else 
     20   USE domqco   , ONLY: dom_qco_zgr      ! vertical scale factor interpolation 
     21#endif 
    1822   USE domngb  , ONLY: dom_ngb          ! find the closest grid point from a given lon/lat position 
    1923   ! 
     
    4347   !! * Substitutions 
    4448#  include "do_loop_substitute.h90" 
     49#  include "domzgr_substitute.h90" 
    4550   !!---------------------------------------------------------------------- 
    4651   !! NEMO/OCE 4.0 , NEMO Consortium (2018) 
     
    112117      vv   (:,:,:,Kbb)   = vv   (:,:,:,Kmm) 
    113118      ssh (:,:,Kbb)     = ssh (:,:,Kmm) 
     119#if ! defined key_qco 
    114120      e3t(:,:,:,Kbb)   = e3t(:,:,:,Kmm) 
    115   
     121#endif  
    116122      ! prepare writing restart 
    117123      IF( lwxios ) THEN 
     
    135141      INTEGER, INTENT(in) :: Kmm    ! ocean time level index 
    136142      !!---------------------------------------------------------------------- 
     143      INTEGER :: jk                               ! loop index 
     144      REAL(wp), DIMENSION(jpi,jpj,jpk) :: ze3t, ze3u, ze3v, zgdepw  ! e3t , e3u, e3v !!st patch to use substitution 
     145      !!---------------------------------------------------------------------- 
     146      ! 
     147      DO jk = 1, jpk 
     148         ze3t(:,:,jk) = e3t(:,:,jk,Kmm) 
     149         ze3u(:,:,jk) = e3u(:,:,jk,Kmm) 
     150         ze3v(:,:,jk) = e3v(:,:,jk,Kmm) 
     151         ! 
     152         zgdepw(:,:,jk) = gdepw(:,:,jk,Kmm) 
     153      END DO  
    137154      ! 
    138155      IF( lwxios ) CALL iom_swap( cwxios_context ) 
    139156      CALL iom_rstput( kt, nitrst, numrow, 'tmask'  , tmask , ldxios = lwxios ) 
    140157      CALL iom_rstput( kt, nitrst, numrow, 'ssmask' , ssmask, ldxios = lwxios ) 
    141       CALL iom_rstput( kt, nitrst, numrow, 'e3t_n'  , e3t(:,:,:,Kmm) , ldxios = lwxios ) 
    142       CALL iom_rstput( kt, nitrst, numrow, 'e3u_n'  , e3u(:,:,:,Kmm) , ldxios = lwxios ) 
    143       CALL iom_rstput( kt, nitrst, numrow, 'e3v_n'  , e3v(:,:,:,Kmm) , ldxios = lwxios ) 
    144       CALL iom_rstput( kt, nitrst, numrow, 'gdepw_n', gdepw(:,:,:,Kmm) , ldxios = lwxios ) 
     158      CALL iom_rstput( kt, nitrst, numrow, 'e3t_n'  , ze3t , ldxios = lwxios ) 
     159      CALL iom_rstput( kt, nitrst, numrow, 'e3u_n'  , ze3u , ldxios = lwxios ) 
     160      CALL iom_rstput( kt, nitrst, numrow, 'e3v_n'  , ze3v , ldxios = lwxios ) 
     161      CALL iom_rstput( kt, nitrst, numrow, 'gdepw_n', zgdepw , ldxios = lwxios ) 
    145162      IF( lwxios ) CALL iom_swap( cxios_context ) 
    146163      ! 
     
    209226      IF(lwp) write(numout,*) 'isfcpl_ssh : recompute scale factor from ssh (new wet cell,Kmm)' 
    210227      IF(lwp) write(numout,*) '~~~~~~~~~~~' 
     228#if ! defined key_qco 
    211229      DO jk = 1, jpk 
    212          e3t(:,:,jk,Kmm) =  e3t_0(:,:,jk) * ( ht_0(:,:) + ssh(:,:,Kmm) ) & 
    213              &                          / ( ht_0(:,:) + 1._wp - ssmask(:,:) ) * tmask(:,:,jk)   & 
    214              &          + e3t_0(:,:,jk)                               * (1._wp -tmask(:,:,jk)) 
     230         e3t(:,:,jk,Kmm) = e3t_0(:,:,jk) * ( 1._wp + (ht_0(:,:) + ssh(:,:,Kmm)) * r1_ht_0(:,:) ) 
    215231      END DO 
    216232      e3t(:,:,:,Kbb) = e3t(:,:,:,Kmm) 
    217233      CALL dom_vvl_zgr(Kbb, Kmm, Kaa) 
     234#else 
     235      CALL dom_qco_zgr(Kbb, Kmm, Kaa) 
     236#endif 
    218237      ! 
    219238   END SUBROUTINE isfcpl_ssh 
     
    400419         ! 1.1: get volume flux before coupling (>0 out) 
    401420         DO_2D_00_00 
    402             zqvolb(ji,jj,jk) =  (   e2u(ji,jj) * ze3u_b(ji,jj,jk) * uu(ji,jj,jk,Kmm) - e2u(ji-1,jj  ) * ze3u_b(ji-1,jj  ,jk) * uu(ji-1,jj  ,jk,Kmm)    & 
    403                &                  + e1v(ji,jj) * ze3v_b(ji,jj,jk) * vv(ji,jj,jk,Kmm) - e1v(ji  ,jj-1) * ze3v_b(ji  ,jj-1,jk) * vv(ji  ,jj-1,jk,Kmm)  ) & 
     421            zqvolb(ji,jj,jk) =    & 
     422               &  (   e2u(ji  ,jj  ) * ze3u_b(ji  ,jj  ,jk) * uu(ji  ,jj  ,jk,Kmm)      & 
     423               &    - e2u(ji-1,jj  ) * ze3u_b(ji-1,jj  ,jk) * uu(ji-1,jj  ,jk,Kmm)      & 
     424               &    + e1v(ji  ,jj  ) * ze3v_b(ji  ,jj  ,jk) * vv(ji  ,jj  ,jk,Kmm)      & 
     425               &    - e1v(ji  ,jj-1) * ze3v_b(ji  ,jj-1,jk) * vv(ji  ,jj-1,jk,Kmm)  )   & 
    404426               &                * ztmask_b(ji,jj,jk) 
    405427         END_2D 
     
    412434         ! compute volume flux divergence after coupling 
    413435         DO_2D_00_00 
    414             zqvoln(ji,jj,jk) = (   e2u(ji,jj) * e3u(ji,jj,jk,Kmm) * uu(ji,jj,jk,Kmm) - e2u(ji-1,jj  ) * e3u(ji-1,jj  ,jk,Kmm) * uu(ji-1,jj  ,jk,Kmm)    & 
    415                &                 + e1v(ji,jj) * e3v(ji,jj,jk,Kmm) * vv(ji,jj,jk,Kmm) - e1v(ji  ,jj-1) * e3v(ji  ,jj-1,jk,Kmm) * vv(ji  ,jj-1,jk,Kmm)  ) & 
     436            zqvoln(ji,jj,jk) =   & 
     437               &  (   e2u(ji  ,jj  ) * e3u(ji  ,jj  ,jk,Kmm) * uu(ji  ,jj  ,jk,Kmm)    & 
     438               &    - e2u(ji-1,jj  ) * e3u(ji-1,jj  ,jk,Kmm) * uu(ji-1,jj  ,jk,Kmm)    & 
     439               &    + e1v(ji  ,jj  ) * e3v(ji  ,jj  ,jk,Kmm) * vv(ji  ,jj  ,jk,Kmm)    & 
     440               &    - e1v(ji  ,jj-1) * e3v(ji  ,jj-1,jk,Kmm) * vv(ji  ,jj-1,jk,Kmm)  ) & 
    416441               &               * tmask(ji,jj,jk) 
    417442         END_2D 
     
    523548 
    524549               ! volume diff 
    525                zdvol = e3t(ji,jj,jk,Kmm) * tmask(ji,jj,jk) - ze3t_b(ji,jj,jk) * ztmask_b(ji,jj,jk) 
     550               zdvol =   e3t  (ji,jj,jk,Kmm) *  tmask  (ji,jj,jk)   & 
     551                  &   - ze3t_b(ji,jj,jk    ) * ztmask_b(ji,jj,jk) 
    526552 
    527553               ! heat diff 
     
    555581            DO ji = nldi,nlei 
    556582               jip1=MIN(ji+1,jpi) ; jim1=MAX(ji-1,1) ; jjp1=MIN(jj+1,jpj) ; jjm1=MAX(jj-1,1) ; 
    557                IF ( tmask(ji,jj,jk) == 0._wp .AND. ztmask_b(ji,jj,jk) == 1._wp ) nisfl(narea) = nisfl(narea) + MAX(SUM(tmask(jim1:jip1,jjm1:jjp1,jk)),1._wp) 
     583               IF ( tmask(ji,jj,jk) == 0._wp .AND. ztmask_b(ji,jj,jk) == 1._wp ) THEN  
     584                  nisfl(narea) = nisfl(narea) + MAX(SUM(tmask(jim1:jip1,jjm1:jjp1,jk)),1._wp) 
     585               ENDIF 
    558586            ENDDO 
    559587         ENDDO 
  • NEMO/trunk/src/OCE/ISF/isfdiags.F90

    r12905 r13237  
    2626   !! * Substitutions 
    2727#  include "do_loop_substitute.h90" 
     28#  include "domzgr_substitute.h90" 
    2829   !!---------------------------------------------------------------------- 
    2930   !! NEMO/OCE 4.0 , NEMO Consortium (2018) 
  • NEMO/trunk/src/OCE/ISF/isfdynatf.F90

    r12489 r13237  
    1414 
    1515   USE phycst , ONLY: r1_rho0         ! physical constant 
    16    USE dom_oce, ONLY: tmask, ssmask, ht, e3t, r1_e1e2t   ! time and space domain 
     16   USE dom_oce                        ! time and space domain 
     17   USE oce, ONLY : ssh                ! sea-surface height !!st needed for substitution 
    1718 
    1819   USE in_out_manager 
     
    2526   !! * Substitutions 
    2627#  include "do_loop_substitute.h90" 
     28#  include "domzgr_substitute.h90" 
    2729 
    2830CONTAINS 
     
    8183      ! add the increment 
    8284      DO jk = 1, jpkm1 
    83          pe3t_f(:,:,jk) = pe3t_f(:,:,jk) - tmask(:,:,jk) * zfwfinc(:,:) * e3t(:,:,jk,Kmm) 
     85         pe3t_f(:,:,jk) = pe3t_f(:,:,jk) - tmask(:,:,jk) * zfwfinc(:,:)   & 
     86            &                              * e3t(:,:,jk,Kmm) 
    8487      END DO 
    8588      ! 
  • NEMO/trunk/src/OCE/ISF/isfhdiv.F90

    r12489 r13237  
    2626   !! * Substitutions 
    2727#  include "do_loop_substitute.h90" 
     28#  include "domzgr_substitute.h90" 
    2829 
    2930CONTAINS 
     
    134135      ! 
    135136      DO jk=1,jpk  
    136          phdiv(:,:,jk) =  phdiv(:,:,jk) + pqvol(:,:,jk) * r1_e1e2t(:,:) / e3t(:,:,jk,Kmm) 
     137         phdiv(:,:,jk) =  phdiv(:,:,jk) + pqvol(:,:,jk) * r1_e1e2t(:,:)   & 
     138            &                             / e3t(:,:,jk,Kmm) 
    137139      END DO 
    138140      ! 
  • NEMO/trunk/src/OCE/ISF/isfload.F90

    r12340 r13237  
    1313   USE isf_oce, ONLY: cn_isfload, rn_isfload_T, rn_isfload_S ! ice shelf variables 
    1414 
    15    USE dom_oce, ONLY: e3w, gdept, risfdep, mikt     ! vertical scale factor 
     15   USE dom_oce                                      ! vertical scale factor 
    1616   USE eosbn2 , ONLY: eos                           ! eos routine 
    1717 
     
    2626   !! * Substitutions 
    2727#  include "do_loop_substitute.h90" 
     28#  include "domzgr_substitute.h90" 
    2829 
    2930CONTAINS 
     
    99100            ! 
    100101            ! top layer of the ice shelf 
    101             pisfload(ji,jj) = pisfload(ji,jj) + (znad + zrhd(ji,jj,1) ) * e3w(ji,jj,1,Kmm) 
     102            pisfload(ji,jj) = pisfload(ji,jj) + (znad + zrhd(ji,jj,1) )   & 
     103               &                                * e3w(ji,jj,1,Kmm) 
    102104            ! 
    103105            ! core layers of the ice shelf 
    104106            DO jk = 2, ikt-1 
    105                pisfload(ji,jj) = pisfload(ji,jj) + (2._wp * znad + zrhd(ji,jj,jk-1) + zrhd(ji,jj,jk)) * e3w(ji,jj,jk,Kmm) 
     107               pisfload(ji,jj) = pisfload(ji,jj) + (2._wp * znad + zrhd(ji,jj,jk-1) + zrhd(ji,jj,jk))   & 
     108                  &                                * e3w(ji,jj,jk,Kmm) 
    106109            END DO 
    107110            ! 
  • NEMO/trunk/src/OCE/ISF/isfstp.F90

    r12242 r13237  
    1313   !!   isfstp       : compute iceshelf melt and heat flux 
    1414   !!---------------------------------------------------------------------- 
    15    ! 
    1615   USE isf_oce                                      ! isf variables 
    1716   USE isfload, ONLY: isf_load                      ! ice shelf load 
     
    2120   USE isfcpl , ONLY: isfcpl_rst_write, isfcpl_init ! isf variables 
    2221 
    23    USE dom_oce, ONLY: ht, e3t, ln_isfcav, ln_linssh     ! ocean space and time domain 
     22   USE dom_oce        ! ocean space and time domain 
     23   USE oce      , ONLY: ssh                           ! sea surface height 
    2424   USE domvvl,  ONLY: ln_vvl_zstar                      ! zstar logical 
    2525   USE zdfdrg,  ONLY: r_Cdmin_top, r_ke0_top            ! vertical physics: top/bottom drag coef. 
     
    3131 
    3232   IMPLICIT NONE 
    33  
    3433   PRIVATE 
    3534 
    3635   PUBLIC   isf_stp, isf_init, isf_nam  ! routine called in sbcmod and divhor 
    3736 
     37   !! * Substitutions 
     38#  include "domzgr_substitute.h90" 
    3839   !!---------------------------------------------------------------------- 
    3940   !! NEMO/OCE 4.0 , NEMO Consortium (2018) 
     
    4142   !! Software governed by the CeCILL license (see ./LICENSE) 
    4243   !!---------------------------------------------------------------------- 
     44 
    4345CONTAINS 
    4446  
     
    6062      INTEGER, INTENT(in) ::   kt   ! ocean time step 
    6163      INTEGER, INTENT(in) ::   Kmm  ! ocean time level index 
     64      !!---------------------------------------------------------------------- 
     65      INTEGER :: jk                               ! loop index 
     66      REAL(wp), DIMENSION(jpi,jpj,jpk) :: ze3t    ! e3t  
    6267      !!--------------------------------------------------------------------- 
    6368      ! 
     
    7883         ! 1.2: compute misfkb, rhisf_tbl, rfrac (deepest level, thickness, fraction of deepest cell affected by tbl) 
    7984         rhisf_tbl_cav(:,:) = rn_htbl * mskisf_cav(:,:) 
    80          CALL isf_tbl_lvl(ht, e3t(:,:,:,Kmm), misfkt_cav, misfkb_cav, rhisf_tbl_cav, rfrac_tbl_cav) 
     85         DO jk = 1, jpk 
     86            ze3t(:,:,jk) = e3t(:,:,jk,Kmm) 
     87         END DO  
     88         CALL isf_tbl_lvl(ht(:,:), ze3t, misfkt_cav, misfkb_cav, rhisf_tbl_cav, rfrac_tbl_cav) 
    8189         ! 
    8290         ! 1.3: compute ice shelf melt 
     
    100108         ! by simplicity, we assume the top level where param applied do not change with time (done in init part) 
    101109         rhisf_tbl_par(:,:) = rhisf0_tbl_par(:,:) 
    102          CALL isf_tbl_lvl(ht, e3t(:,:,:,Kmm), misfkt_par, misfkb_par, rhisf_tbl_par, rfrac_tbl_par) 
     110         DO jk = 1, jpk 
     111            ze3t(:,:,jk) = e3t(:,:,jk,Kmm) 
     112         END DO 
     113         CALL isf_tbl_lvl(ht(:,:), ze3t, misfkt_par, misfkb_par, rhisf_tbl_par, rfrac_tbl_par) 
    103114         ! 
    104115         ! 2.3: compute ice shelf melt 
  • NEMO/trunk/src/OCE/ISF/isftbl.F90

    r12340 r13237  
    2525   !! * Substitutions 
    2626#  include "do_loop_substitute.h90" 
     27#  include "domzgr_substitute.h90" 
    2728 
    2829CONTAINS 
     
    5657      REAL(wp), DIMENSION(jpi,jpj) :: zhtbl   ! thickness of the tbl 
    5758      REAL(wp), DIMENSION(jpi,jpj) :: zfrac   ! thickness of the tbl 
     59      INTEGER :: jk                            ! loop index 
     60      REAL(wp), DIMENSION(jpi,jpj,jpk) :: ze3t,ze3u,ze3v ! e3  
    5861      !!-------------------------------------------------------------------- 
    5962      !  
     
    6467         zhtbl = phtbl 
    6568         ! 
     69         DO jk = 1, jpk 
     70            ze3u(:,:,jk) = e3u(:,:,jk,Kmm) 
     71         END DO  
    6672         ! compute tbl lvl and thickness 
    67          CALL isf_tbl_lvl( hu(:,:,Kmm), e3u(:,:,:,Kmm), ktop, ikbot, zhtbl, zfrac ) 
     73         CALL isf_tbl_lvl( hu(:,:,Kmm), ze3u, ktop, ikbot, zhtbl, zfrac ) 
    6874         ! 
    6975         ! compute tbl property at U point 
    70          CALL isf_tbl_avg( miku, ikbot, zhtbl, zfrac, e3u(:,:,:,Kmm), pvarin, zvarout ) 
     76         CALL isf_tbl_avg( miku, ikbot, zhtbl, zfrac, ze3u, pvarin, zvarout ) 
    7177         ! 
    7278         ! compute tbl property at T point 
     
    8288         zhtbl = phtbl 
    8389         ! 
     90         DO jk = 1, jpk 
     91            ze3v(:,:,jk) = e3v(:,:,jk,Kmm) 
     92         END DO  
    8493         ! compute tbl lvl and thickness 
    85          CALL isf_tbl_lvl( hv(:,:,Kmm), e3v(:,:,:,Kmm), ktop, ikbot, zhtbl, zfrac ) 
     94         CALL isf_tbl_lvl( hv(:,:,Kmm), ze3v, ktop, ikbot, zhtbl, zfrac ) 
    8695         ! 
    8796         ! compute tbl property at V point 
    88          CALL isf_tbl_avg( mikv, ikbot, zhtbl, zfrac, e3v(:,:,:,Kmm), pvarin, zvarout ) 
     97         CALL isf_tbl_avg( mikv, ikbot, zhtbl, zfrac, ze3v, pvarin, zvarout ) 
    8998         ! 
    9099         ! pvarout is an averaging of wet point 
     
    98107         ! 
    99108         ! compute tbl property at T point 
    100          CALL isf_tbl_avg( ktop, kbot, phtbl, pfrac, e3t(:,:,:,Kmm), pvarin, pvarout ) 
     109         DO jk = 1, jpk 
     110            ze3t(:,:,jk) = e3t(:,:,jk,Kmm) 
     111         END DO  
     112         CALL isf_tbl_avg( ktop, kbot, phtbl, pfrac, ze3t, pvarin, pvarout ) 
    101113         ! 
    102114      END SELECT 
     
    212224      ! phtbl need to be bounded by water column thickness before 
    213225      ! test: if htbl = water column thickness, should return mbathy 
    214       ! test: if htbl = 0 should return ktop (phtbl cap to e3t(ji,jj,1)) 
     226      ! test: if htbl = 0 should return ktop (phtbl cap to pe3t(ji,jj,1)) 
    215227      ! 
    216228      ! get ktbl 
Note: See TracChangeset for help on using the changeset viewer.