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 8930 for branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC – NEMO

Ignore:
Timestamp:
2017-12-06T18:17:10+01:00 (6 years ago)
Author:
acc
Message:

Branch 2017/dev_CNRS_2017. Merge in NOC-originated OSMOSIS additions from dev_r7881_ENHANCE09_RK3. SETTE tested successfully apart from WISOMIP_ST (expected) and WORCA2AGR_ST (which differs only after 115 timesteps and seems to be a local issue)

Location:
branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC
Files:
1 added
10 edited

Legend:

Unmodified
Added
Removed
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC/C1D/dyncor_c1d.F90

    r7646 r8930  
    2121   USE prtctl         ! Print control 
    2222 
     23   USE sbcwave        ! Surface Waves (add Stokes-Coriolis force) 
     24   USE sbc_oce , ONLY : ln_stcor    ! use Stoke-Coriolis force 
     25    
    2326   IMPLICIT NONE 
    2427   PRIVATE 
     
    7174      ENDIF 
    7275      ! 
    73       DO jk = 1, jpkm1 
    74          DO jj = 2, jpjm1 
    75             DO ji = fs_2, fs_jpim1   ! vector opt. 
    76                ua(ji,jj,jk) = ua(ji,jj,jk) + ff_t(ji,jj) * vn(ji,jj,jk) 
    77                va(ji,jj,jk) = va(ji,jj,jk) - ff_t(ji,jj) * un(ji,jj,jk) 
     76      IF( ln_stcor ) THEN 
     77         DO jk = 1, jpkm1 
     78            DO jj = 2, jpjm1 
     79               DO ji = fs_2, fs_jpim1   ! vector opt. 
     80                  ua(ji,jj,jk) = ua(ji,jj,jk) + ff_t(ji,jj) * (vn(ji,jj,jk) + vsd(ji,jj,jk)) 
     81                  va(ji,jj,jk) = va(ji,jj,jk) - ff_t(ji,jj) * (un(ji,jj,jk) + usd(ji,jj,jk)) 
     82               END DO 
    7883            END DO 
    7984         END DO 
    80       END DO    
     85      ELSE 
     86         DO jk = 1, jpkm1 
     87            DO jj = 2, jpjm1 
     88               DO ji = fs_2, fs_jpim1   ! vector opt. 
     89                  ua(ji,jj,jk) = ua(ji,jj,jk) + ff_t(ji,jj) * vn(ji,jj,jk) 
     90                  va(ji,jj,jk) = va(ji,jj,jk) - ff_t(ji,jj) * un(ji,jj,jk) 
     91               END DO 
     92            END DO 
     93         END DO 
     94      END IF 
     95       
    8196      ! 
    8297      IF(ln_ctl)   CALL prt_ctl( tab3d_1=ua, clinfo1=' cor  - Ua: ', mask1=umask,  & 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC/C1D/step_c1d.F90

    r8882 r8930  
    7878                         CALL zdf_phy( kstp )         ! vertical physics update (bfr, avt, avs, avm + MLD) 
    7979 
     80      IF(.NOT.ln_linssh )   CALL ssh_nxt       ( kstp )  ! after ssh (includes call to div_hor) 
     81      IF(.NOT.ln_linssh )   CALL dom_vvl_sf_nxt( kstp )  ! after vertical scale factors  
     82 
     83      IF(.NOT.ln_linssh )   CALL wzv           ( kstp )  ! now cross-level velocity  
    8084      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
    8185      ! diagnostics and outputs             (ua, va, ta, sa used as workspace) 
     
    100104      IF( ln_traqsr )   CALL tra_qsr( kstp )       ! penetrative solar radiation qsr 
    101105      IF( ln_tradmp )   CALL tra_dmp( kstp )       ! internal damping trends- tracers 
     106      IF(.NOT.ln_linssh)CALL tra_adv( kstp )       ! horizontal & vertical advection 
     107      IF( lk_zdfosm  )  CALL tra_osm( kstp )       ! OSMOSIS non-local tracer fluxes 
    102108                        CALL tra_zdf( kstp )       ! vertical mixing 
    103109                        CALL eos( tsn, rhd, rhop, gdept_0(:,:,:) )   ! now potential density for zdfmxl 
     
    115121      IF( ln_dyndmp )   CALL dyn_dmp    ( kstp )   ! internal damping trends- momentum 
    116122                        CALL dyn_cor_c1d( kstp )   ! vorticity term including Coriolis 
     123      IF( lk_zdfosm  )  CALL dyn_osm    ( kstp )   ! OSMOSIS non-local velocity fluxes 
    117124                        CALL dyn_zdf    ( kstp )   ! vertical diffusion 
    118125                        CALL dyn_nxt    ( kstp )   ! lateral velocity at next time step 
     126      IF(.NOT.ln_linssh)CALL ssh_swp    ( kstp )   ! swap of sea surface height 
     127 
     128      IF(.NOT.ln_linssh)CALL dom_vvl_sf_swp( kstp )! swap of vertical scale factors 
    119129 
    120130      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC/DIA/diacfl.F90

    r8882 r8930  
    105105      !                    ! write out to file 
    106106      IF( lwp ) THEN 
    107          WRITE(numcfl,FMT='(2x,i4,5x,a6,5x,f6.4,1x,i4,1x,i4,1x,i4)') kt, 'Max Cu', zCu_max, iloc_u(1), iloc_u(2), iloc_u(3) 
    108          WRITE(numcfl,FMT='(11x,     a6,5x,f6.4,1x,i4,1x,i4,1x,i4)')     'Max Cv', zCv_max, iloc_v(1), iloc_v(2), iloc_v(3) 
    109          WRITE(numcfl,FMT='(11x,     a6,5x,f6.4,1x,i4,1x,i4,1x,i4)')     'Max Cw', zCw_max, iloc_w(1), iloc_w(2), iloc_w(3) 
     107         WRITE(numcfl,FMT='(2x,i4,5x,a6,4x,f7.4,1x,i4,1x,i4,1x,i4)') kt, 'Max Cu', zCu_max, iloc_u(1), iloc_u(2), iloc_u(3) 
     108         WRITE(numcfl,FMT='(11x,     a6,4x,f7.4,1x,i4,1x,i4,1x,i4)')     'Max Cv', zCv_max, iloc_v(1), iloc_v(2), iloc_v(3) 
     109         WRITE(numcfl,FMT='(11x,     a6,4x,f7.4,1x,i4,1x,i4,1x,i4)')     'Max Cw', zCw_max, iloc_w(1), iloc_w(2), iloc_w(3) 
    110110      ENDIF 
    111111      ! 
     
    119119         ! to ascii file 
    120120         WRITE(numcfl,*) '******************************************' 
    121          WRITE(numcfl,FMT='(3x,a12,7x,f6.4,1x,i4,1x,i4,1x,i4)') 'Run Max Cu', rCu_max, nCu_loc(1), nCu_loc(2), nCu_loc(3) 
     121         WRITE(numcfl,FMT='(3x,a12,6x,f7.4,1x,i4,1x,i4,1x,i4)') 'Run Max Cu', rCu_max, nCu_loc(1), nCu_loc(2), nCu_loc(3) 
    122122         WRITE(numcfl,FMT='(3x,a8,11x,f15.1)') ' => dt/C', z2dt/rCu_max 
    123123         WRITE(numcfl,*) '******************************************' 
    124          WRITE(numcfl,FMT='(3x,a12,7x,f6.4,1x,i4,1x,i4,1x,i4)') 'Run Max Cv', rCv_max, nCv_loc(1), nCv_loc(2), nCv_loc(3) 
     124         WRITE(numcfl,FMT='(3x,a12,6x,f7.4,1x,i4,1x,i4,1x,i4)') 'Run Max Cv', rCv_max, nCv_loc(1), nCv_loc(2), nCv_loc(3) 
    125125         WRITE(numcfl,FMT='(3x,a8,11x,f15.1)') ' => dt/C', z2dt/rCv_max 
    126126         WRITE(numcfl,*) '******************************************' 
    127          WRITE(numcfl,FMT='(3x,a12,7x,f6.4,1x,i4,1x,i4,1x,i4)') 'Run Max Cw', rCw_max, nCw_loc(1), nCw_loc(2), nCw_loc(3) 
     127         WRITE(numcfl,FMT='(3x,a12,6x,f7.4,1x,i4,1x,i4,1x,i4)') 'Run Max Cw', rCw_max, nCw_loc(1), nCw_loc(2), nCw_loc(3) 
    128128         WRITE(numcfl,FMT='(3x,a8,11x,f15.1)') ' => dt/C', z2dt/rCw_max 
    129129         CLOSE( numcfl )  
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC/DIA/diahth.F90

    r6140 r8930  
    139139      DO jj = 1, jpj 
    140140         DO ji = 1, jpi 
    141             zztmp = bathy(ji,jj) 
     141            zztmp = gdepw_n(ji,jj,mbkt(ji,jj)+1)  
    142142            hth     (ji,jj) = zztmp 
    143143            zabs2   (ji,jj) = zztmp 
     
    150150         DO jj = 1, jpj 
    151151            DO ji = 1, jpi 
    152                zztmp = bathy(ji,jj) 
     152               zztmp = gdepw_n(ji,jj,mbkt(ji,jj)+1)  
    153153               zrho0_3(ji,jj) = zztmp 
    154154               zrho0_1(ji,jj) = zztmp 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC/DOM/domain.F90

    r8882 r8930  
    3737   USE domvvl         ! variable volume 
    3838   USE c1d            ! 1D configuration 
    39    USE domc1d         ! 1D configuration: column location 
    4039   USE dyncor_c1d     ! 1D configuration: Coriolis term    (cor_c1d routine) 
    4140   USE wet_dry        ! wetting and drying 
     
    117116         WRITE(numout,*)     '              nn_cfg = ', nn_cfg 
    118117      ENDIF 
    119       ! 
    120       !       
    121 !!gm  This should be removed with the new configuration interface 
    122       IF( lk_c1d .AND. ln_c1d_locpt )  CALL dom_c1d( rn_lat1d, rn_lon1d ) 
    123 !!gm end 
    124118      ! 
    125119      !           !==  Reference coordinate system  ==! 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC/SBC/sbcwave.F90

    r8882 r8930  
    8787      INTEGER  ::   jj, ji, jk   ! dummy loop argument 
    8888      INTEGER  ::   ik           ! local integer  
    89       REAL(wp) ::  ztransp, zfac, ztemp, zsp0 
    90       REAL(wp) ::  zdep_u, zdep_v, zkh_u, zkh_v, zda_u, zda_v 
     89      REAL(wp) ::  ztransp, zfac, zsp0 
     90      REAL(wp) ::  zdepth, zsqrt_depth,  zexp_depth, z_two_thirds, zsqrtpi !sqrt of pi 
     91      REAL(wp) ::  zbot_u, zbot_v, zkb_u, zkb_v, zke3_u, zke3_v, zda_u, zda_v 
     92      REAL(wp) ::  zstokes_psi_u_bot, zstokes_psi_v_bot 
    9193      REAL(wp), DIMENSION(:,:)  , POINTER ::   zk_t, zk_u, zk_v, zu0_sd, zv0_sd   ! 2D workspace 
     94      REAL(wp), DIMENSION(:,:)  , POINTER ::   zstokes_psi_u_top, zstokes_psi_v_top   ! 2D workspace 
    9295      REAL(wp), DIMENSION(:,:,:), POINTER ::   ze3divh                            ! 3D workspace 
    9396      !!--------------------------------------------------------------------- 
     
    9598      CALL wrk_alloc( jpi,jpj,jpk,   ze3divh ) 
    9699      CALL wrk_alloc( jpi,jpj,       zk_t, zk_u, zk_v, zu0_sd, zv0_sd ) 
    97       ! 
    98       ! 
     100      CALL wrk_alloc( jpi,jpj,       zstokes_psi_u_top, zstokes_psi_v_top) 
     101      ! 
     102      ! 
     103      zsqrtpi = SQRT(rpi) 
     104      z_two_thirds = 2.0_wp / 3.0_wp 
    99105      zfac =  2.0_wp * rpi / 16.0_wp 
    100106      DO jj = 1, jpj                ! exp. wave number at t-point    (Eq. (19) in Breivick et al. (2014) ) 
     
    106112               tsd2d(ji,jj) = zsp0 
    107113               ! Wavenumber scale 
    108                zk_t(ji,jj) = ABS( zsp0 ) / MAX( ABS( 5.97_wp*ztransp ) , 0.0000001_wp ) 
    109          END DO 
    110       END DO       
     114               zk_t(ji,jj) =  (1.0_wp-2.0_wp/3.0_wp)*zsp0/MAX(2.0_wp*ztransp,0.0000001_wp) 
     115         END DO 
     116      END DO 
     117      ! 
    111118      DO jj = 1, jpjm1              ! exp. wave number & Stokes drift velocity at u- & v-points 
    112119         DO ji = 1, jpim1 
     
    120127      ! 
    121128      !                       !==  horizontal Stokes Drift 3D velocity  ==! 
     129      DO jj = 1, jpjm1              ! exp. wave number & Stokes drift velocity at u- & v-points 
     130         DO ji = 1, jpim1 
     131            zstokes_psi_u_top(ji,jj) = 0._wp 
     132            zstokes_psi_v_top(ji,jj) = 0._wp 
     133         END DO 
     134      END DO 
     135 
    122136      DO jk = 1, jpkm1 
    123137         DO jj = 2, jpjm1 
    124138            DO ji = 2, jpim1 
    125                zdep_u = 0.5_wp * ( gdept_n(ji,jj,jk) + gdept_n(ji+1,jj,jk) ) 
    126                zdep_v = 0.5_wp * ( gdept_n(ji,jj,jk) + gdept_n(ji,jj+1,jk) ) 
    127                !                           
    128                zkh_u = zk_u(ji,jj) * zdep_u     ! k * depth 
    129                zkh_v = zk_v(ji,jj) * zdep_v 
    130                !                                ! Depth attenuation 
    131                zda_u = EXP( -2.0_wp*zkh_u ) / ( 1.0_wp + 8.0_wp*zkh_u ) 
    132                zda_v = EXP( -2.0_wp*zkh_v ) / ( 1.0_wp + 8.0_wp*zkh_v ) 
     139               zbot_u = 0.5_wp * ( gdepw_n(ji,jj,jk+1) + gdepw_n(ji+1,jj,jk+1) ) 
     140               zbot_v = 0.5_wp * ( gdepw_n(ji,jj,jk+1) + gdepw_n(ji,jj+1,jk+1) ) 
     141               zkb_u = 2.0_wp * zk_u(ji,jj) * zbot_u     ! 2k * bottom depth 
     142               zkb_v = 2.0_wp * zk_v(ji,jj) * zbot_v     ! 2k * bottom depth 
    133143               ! 
     144               zke3_u = MAX(1.e-8_wp, 2.0_wp * zk_u(ji,jj) * e3u_n(ji,jj,jk))     ! 2k * thickness 
     145               zke3_v = MAX(1.e-8_wp, 2.0_wp * zk_v(ji,jj) * e3v_n(ji,jj,jk))     ! 2k * thickness 
     146 
     147               ! Depth attenuation .... do u component first.. 
     148               zdepth=zkb_u 
     149               zsqrt_depth = SQRT(zdepth) 
     150               zexp_depth = EXP(-zdepth) 
     151               zstokes_psi_u_bot = 1.0_wp - zexp_depth  & 
     152                    & - z_two_thirds*(zsqrtpi*zsqrt_depth*zdepth*ERFC(zsqrt_depth) & 
     153                    &                  + 1.0_wp - (1.0_wp + zdepth)*zexp_depth) 
     154               zda_u = (zstokes_psi_u_bot - zstokes_psi_u_top(ji,jj))/zke3_u 
     155               zstokes_psi_u_top(ji,jj) = zstokes_psi_u_bot 
     156 
     157               !         ... and then v component 
     158               zdepth=zkb_v 
     159               zsqrt_depth = SQRT(zdepth) 
     160               zexp_depth = EXP(-zdepth) 
     161               zstokes_psi_v_bot = 1.0_wp - zexp_depth  & 
     162                    & - z_two_thirds*(zsqrtpi*zsqrt_depth*zdepth*ERFC(zsqrt_depth) & 
     163                    &                  + 1.0_wp - (1.0_wp + zdepth)*zexp_depth) 
     164               zda_v = (zstokes_psi_v_bot - zstokes_psi_v_top(ji,jj))/zke3_v 
     165               zstokes_psi_v_top(ji,jj) = zstokes_psi_v_bot 
     166              ! 
    134167               usd(ji,jj,jk) = zda_u * zu0_sd(ji,jj) * umask(ji,jj,jk) 
    135168               vsd(ji,jj,jk) = zda_v * zv0_sd(ji,jj) * vmask(ji,jj,jk) 
     
    190223      CALL wrk_dealloc( jpi,jpj,jpk,   ze3divh ) 
    191224      CALL wrk_dealloc( jpi,jpj,       zk_t, zk_u, zk_v, zu0_sd, zv0_sd ) 
     225      CALL wrk_dealloc( jpi,jpj,       zstokes_psi_u_top, zstokes_psi_v_top) 
    192226      ! 
    193227   END SUBROUTINE sbc_stokes 
     
    347381         ALLOCATE( div_sd(jpi,jpj) ) 
    348382         ALLOCATE( tsd2d (jpi,jpj) ) 
     383 
     384         ut0sd(:,:) = 0._wp 
     385         vt0sd(:,:) = 0._wp 
     386         hsw(:,:) = 0._wp 
     387         wmp(:,:) = 0._wp 
     388 
    349389         usd(:,:,:) = 0._wp 
    350390         vsd(:,:,:) = 0._wp 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC/ZDF/zdf_oce.F90

    r8882 r8930  
    2222   LOGICAL , PUBLIC ::   ln_zdfric   !: Richardson depend coefficients 
    2323   LOGICAL , PUBLIC ::   ln_zdftke   !: Turbulent Kinetic Energy closure 
    24    LOGICAL , PUBLIC ::   ln_zdfgls   !: Generic Length Sclare closure 
     24   LOGICAL , PUBLIC ::   ln_zdfgls   !: Generic Length Scale closure 
     25   LOGICAL , PUBLIC ::   ln_zdfosm   !: OSMOSIS BL closure 
    2526   !                             ! convection 
    2627   LOGICAL , PUBLIC ::   ln_zdfevd   !: convection: enhanced vertical diffusion flag 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC/ZDF/zdfphy.F90

    r8882 r8930  
    1818   USE zdftke         ! vertical physics: TKE vertical mixing 
    1919   USE zdfgls         ! vertical physics: GLS vertical mixing 
     20   USE zdfosm         ! vertical physics: OSMOSIS vertical mixing 
    2021   USE zdfddm         ! vertical physics: double diffusion mixing       
    2122   USE zdfevd         ! vertical physics: convection via enhanced vertical diffusion   
     
    4950   INTEGER, PARAMETER ::   np_TKE = 3   ! Turbulente Kinetic Eenergy closure scheme for Kz 
    5051   INTEGER, PARAMETER ::   np_GLS = 4   ! Generic Length Scale closure scheme for Kz 
     52   INTEGER, PARAMETER ::   np_OSM = 5   ! OSMOSIS-OBL closure scheme for Kz 
    5153 
    5254   LOGICAL ::   l_zdfsh2   ! shear production term flag (=F for CST, =T otherwise (i.e. TKE, GLS, RIC)) 
     
    7274      !! 
    7375      NAMELIST/namzdf/ ln_zdfcst, ln_zdfric, ln_zdftke, ln_zdfgls,   &     ! type of closure scheme 
     76         &             ln_zdfosm,                                    &     ! type of closure scheme 
    7477         &             ln_zdfevd, nn_evdm, rn_evd ,                  &     ! convection : evd 
    7578         &             ln_zdfnpc, nn_npc , nn_npcp,                  &     ! convection : npc 
     
    102105         WRITE(numout,*) '         Turbulent Kinetic Energy closure (TKE)  ln_zdftke = ', ln_zdftke 
    103106         WRITE(numout,*) '         Generic Length Scale closure (GLS)      ln_zdfgls = ', ln_zdfgls 
     107         WRITE(numout,*) '         OSMOSIS-OBL closure (OSM)               ln_zdfosm = ', ln_zdfosm 
    104108         WRITE(numout,*) '      convection: ' 
    105109         WRITE(numout,*) '         enhanced vertical diffusion             ln_zdfevd = ', ln_zdfevd 
     
    156160      ! 
    157161      IF( ln_zdfnpc .AND. ln_zdfevd )   CALL ctl_stop( 'zdf_phy_init: chose between ln_zdfnpc and ln_zdfevd' ) 
     162      IF( ln_zdfosm .AND. ln_zdfevd )   CALL ctl_stop( 'zdf_phy_init: chose between ln_zdfosm and ln_zdfevd' ) 
    158163      IF( lk_top    .AND. ln_zdfnpc )   CALL ctl_stop( 'zdf_phy_init: npc scheme is not working with key_top' ) 
     164      IF( lk_top    .AND. ln_zdfosm )   CALL ctl_stop( 'zdf_phy_init: osmosis scheme is not working with key_top' ) 
    159165      IF(lwp) THEN 
    160166         WRITE(numout,*) 
     
    178184      IF( ln_zdftke ) THEN   ;   ioptio = ioptio + 1   ;    nzdf_phy = np_TKE   ;   CALL zdf_tke_init   ;   ENDIF 
    179185      IF( ln_zdfgls ) THEN   ;   ioptio = ioptio + 1   ;    nzdf_phy = np_GLS   ;   CALL zdf_gls_init   ;   ENDIF 
     186      IF( ln_zdfosm ) THEN   ;   ioptio = ioptio + 1   ;    nzdf_phy = np_OSM   ;   CALL zdf_osm_init   ;   ENDIF 
    180187      ! 
    181188      IF( ioptio /= 1 )    CALL ctl_stop( 'zdf_phy_init: one and only one vertical diffusion option has to be defined ' ) 
     
    247254      CASE( np_TKE )   ;   CALL zdf_tke( kt         , zsh2, avm_k, avt_k )    ! TKE closure scheme for Kz 
    248255      CASE( np_GLS )   ;   CALL zdf_gls( kt         , zsh2, avm_k, avt_k )    ! GLS closure scheme for Kz 
     256      CASE( np_OSM )   ;   CALL zdf_osm( kt               , avm_k, avt_k )    ! OSMOSIS closure scheme for Kz 
    249257!     CASE( np_CST )                                  ! Constant Kz (reset avt, avm to the background value) 
    250258!         ! avt_k and avm_k set one for all at initialisation phase 
     
    290298      CALL lbc_lnk( avm  , 'W', 1. )                  ! needed to compute avm at u- and v-points 
    291299      CALL lbc_lnk( avt  , 'W', 1. )                  !!gm  a priori only avm_k and avm are required 
    292       CALL lbc_lnk( avs  , 'W', 1. )                  !!gm  for calculation, keeped here for output only 
     300      CALL lbc_lnk( avs  , 'W', 1. )                  !!gm  for calculation, kept here for output only 
    293301      ! 
    294302      IF( l_zdfdrg ) THEN     ! drag  have been updated (non-linear cases) 
     
    303311         IF( ln_zdfgls )   CALL gls_rst( kt, 'WRITE' ) 
    304312         IF( ln_zdfric )   CALL ric_rst( kt, 'WRITE' )  
     313         ! NB. OSMOSIS restart (osm_rst) will be called in step.F90 after wn has been updated 
    305314      ENDIF 
    306315      ! 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC/step.F90

    r8882 r8930  
    186186                         CALL dyn_vor       ( kstp )  ! vorticity term including Coriolis 
    187187                         CALL dyn_ldf       ( kstp )  ! lateral mixing 
     188      IF( ln_zdfosm  )   CALL dyn_osm       ( kstp )  ! OSMOSIS non-local velocity fluxes 
    188189                         CALL dyn_hpg       ( kstp )  ! horizontal gradient of Hydrostatic pressure 
    189190                         CALL dyn_spg       ( kstp )  ! surface pressure gradient 
     
    243244#endif 
    244245                         CALL tra_adv       ( kstp )  ! horizontal & vertical advection 
     246      IF( ln_zdfosm  )   CALL tra_osm       ( kstp )  ! OSMOSIS non-local tracer fluxes 
     247      IF( lrst_oce .AND. ln_zdfosm ) & 
     248           &             CALL osm_rst( kstp, 'WRITE' )! write OSMOSIS outputs + wn (so must do here) to restarts 
    245249                         CALL tra_ldf       ( kstp )  ! lateral mixing 
    246250 
  • branches/2017/dev_CNRS_2017/NEMOGCM/NEMO/OPA_SRC/step_oce.F90

    r8882 r8930  
    6565 
    6666   USE zdfphy          ! vertical physics manager      (zdf_phy_init routine) 
     67   USE zdfosm  , ONLY : osm_rst, dyn_osm, tra_osm      ! OSMOSIS routines used in step.F90 
    6768 
    6869   USE step_diu        ! Time stepping for diurnal sst 
Note: See TracChangeset for help on using the changeset viewer.