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 11180 for NEMO/branches – NEMO

Changeset 11180 for NEMO/branches


Ignore:
Timestamp:
2019-06-25T18:50:27+02:00 (5 years ago)
Author:
clne
Message:

Initial commit of code for 2d (surge) work in NEMO4.
This is aiming to replicate the 3.6 version in branches/UKMO/dev_r5518_Surge_Modelling

Location:
NEMO/branches/UKMO/NEMO_4.0_surge
Files:
31 added
18 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/UKMO/NEMO_4.0_surge/cfgs/SHARED/field_def_nemo-oce.xml

    r10823 r11180  
    259259          <field id="taum"         long_name="wind stress module"                    standard_name="magnitude_of_surface_downward_stress"                                 unit="N/m2"                           /> 
    260260          <field id="wspd"         long_name="wind speed module"                     standard_name="wind_speed"                                                           unit="m/s"                            /> 
     261          <field id="uwnd"         long_name="u component of wind"       unit="m/s"          grid_ref="grid_U_2D"     /> 
     262          <field id="vwnd"         long_name="v component of wind"       unit="m/s"          grid_ref="grid_V_2D"    /> 
    261263           
    262264          <!-- * variable relative to atmospheric pressure forcing : available with ln_apr_dyn --> 
  • NEMO/branches/UKMO/NEMO_4.0_surge/cfgs/SHARED/namelist_ref

    r10888 r11180  
    7474   ! 
    7575   ln_crs      = .false.   !  Logical switch for coarsening module      (T => fill namcrs) 
     76   ! 
     77   ln_2d        = .false.  !  (=T) run in 2D barotropic mode (no tracer processes or vertical diffusion)  
    7678   ! 
    7779   ln_meshmask = .false.   !  =T create a mesh file 
  • NEMO/branches/UKMO/NEMO_4.0_surge/src/OCE/DIA/diawri.F90

    r10888 r11180  
    185185         CALL iom_put( "taubot", z2d )            
    186186      ENDIF 
    187           
     187  
     188      IF( iom_use("uwnd") ) CALL iom_put( "uwnd" ,   uwnd  )  
     189      IF( iom_use("vwnd") ) CALL iom_put( "vwnd" ,   vwnd  )  
     190           
    188191      CALL iom_put( "uoce", un(:,:,:) )            ! 3D i-current 
    189192      CALL iom_put(  "ssu", un(:,:,1) )            ! surface i-current 
  • NEMO/branches/UKMO/NEMO_4.0_surge/src/OCE/DOM/dom_oce.F90

    r10888 r11180  
    3838   LOGICAL , PUBLIC ::   ln_iscpl       !: coupling with ice sheet 
    3939   LOGICAL , PUBLIC ::   ln_crs         !: Apply grid coarsening to dynamical model output or online passive tracers 
    40  
     40   LOGICAL , PUBLIC ::   ln_2d          !  Default False. If True run in 2D barotropic mode (no tracer processes or vertical diffusion)  
     41  
    4142   !! Free surface parameters 
    4243   !! ======================= 
  • NEMO/branches/UKMO/NEMO_4.0_surge/src/OCE/DOM/domain.F90

    r10888 r11180  
    293293         &             nn_stock, nn_write , ln_mskland  , ln_clobber   , nn_chunksz, nn_euler  ,     & 
    294294         &             ln_cfmeta, ln_iscpl, ln_xios_read, nn_wxios 
    295       NAMELIST/namdom/ ln_linssh, rn_isfhmin, rn_rdt, rn_atfp, ln_crs, ln_meshmask 
     295      NAMELIST/namdom/ ln_linssh, rn_isfhmin, rn_rdt, rn_atfp, ln_crs, ln_meshmask, ln_2d 
    296296#if defined key_netcdf4 
    297297      NAMELIST/namnc4/ nn_nchunks_i, nn_nchunks_j, nn_nchunks_k, ln_nc4zip 
     
    416416         WRITE(numout,*) '      asselin time filter parameter           rn_atfp     = ', rn_atfp 
    417417         WRITE(numout,*) '      online coarsening of dynamical fields   ln_crs      = ', ln_crs 
     418         WRITE(numout,*) '      2D mode                                 ln_2d      = ', ln_2d  
     419         IF(ln_2d) WRITE(numout,*) '   2D mode active: All tracer processes and vertical diffusion turned off'       
    418420      ENDIF 
    419421      ! 
  • NEMO/branches/UKMO/NEMO_4.0_surge/src/OCE/DOM/domzgr.F90

    r10888 r11180  
    134134      IF( ln_sco      )   ioptio = ioptio + 1 
    135135      IF( ioptio /= 1 )   CALL ctl_stop( ' none or several vertical coordinate options used' ) 
    136  
     136      IF (ln_2d .AND. .NOT.ln_sco)   CALL ctl_stop( ' 2D mode must be used with ln_sco' )  
    137137 
    138138      !                                ! top/bottom ocean level indices for t-, u- and v-points (f-point also for top) 
  • NEMO/branches/UKMO/NEMO_4.0_surge/src/OCE/DOM/dtatsd.F90

    r10888 r11180  
    8787            WRITE(numout,*) '   ===>>   T & S data not used' 
    8888         ENDIF 
     89         IF( ln_2d ) WRITE(numout,*) '   2D ocean - ocean will be started at rest and T&S = arbitrary constants'  
    8990      ENDIF 
    9091      ! 
  • NEMO/branches/UKMO/NEMO_4.0_surge/src/OCE/DOM/istate.F90

    r10888 r11180  
    9696         !                                    ! Initialization of ocean to zero 
    9797         ! 
    98          IF( ln_tsd_init ) THEN                
     98         IF (ln_2d) THEN  
     99            IF(lwp) WRITE(numout,*) 'istate_init : 2D case, setting tracers to contants and ocean at rest'  
     100            tsb(:,:,:,:)= 15._wp               ! No tracers, but can't set salinity to 0 otherwise it triggers crash message 
     101            sshb(:,:)   = 0._wp                ! set the ocean at rest 
     102            ub  (:,:,:) = 0._wp                 
     103            vb  (:,:,:) = 0._wp            
     104         ELSE IF( ln_tsd_init ) THEN    
    99105            CALL dta_tsd( nit000, tsb )       ! read 3D T and S data at nit000 
    100106            ! 
  • NEMO/branches/UKMO/NEMO_4.0_surge/src/OCE/DYN/dynhpg.F90

    r10888 r11180  
    493493      END IF 
    494494 
     495    IF (ln_2d) THEN  
     496        ! Surface value  
     497        DO jj = 2, jpjm1  
     498           DO ji = fs_2, fs_jpim1   ! vector opt.  
     499              ! hydrostatic pressure gradient along s-surfaces  
     500              zhpi(ji,jj,1) = zcoef0 * (  e3w_n(ji+1,jj  ,1) * znad   &  
     501                 &                      - e3w_n(ji  ,jj  ,1) * znad ) * r1_e1u(ji,jj)  
     502              zhpj(ji,jj,1) = zcoef0 * (  e3w_n(ji  ,jj+1,1) * znad   &  
     503                 &                      - e3w_n(ji  ,jj  ,1) * znad ) * r1_e2v(ji,jj)  
     504              ! s-coordinate pressure gradient correction  
     505              zuap = -zcoef0 * ( 2._wp * znad )   &  
     506                 &           * ( gde3w_n(ji+1,jj,1) - gde3w_n(ji,jj,1) ) * r1_e1u(ji,jj)  
     507              zvap = -zcoef0 * ( 2._wp * znad )   &  
     508                 &           * ( gde3w_n(ji,jj+1,1) - gde3w_n(ji,jj,1) ) * r1_e2v(ji,jj)  
     509  
     510  
     511              IF( ln_wd_il ) THEN  
     512  
     513                zhpi(ji,jj,1) = zhpi(ji,jj,1) * zcpx(ji,jj)  
     514                zhpj(ji,jj,1) = zhpj(ji,jj,1) * zcpy(ji,jj)   
     515                zuap = zuap * zcpx(ji,jj)  
     516                zvap = zvap * zcpy(ji,jj)  
     517              ENDIF  
     518  
     519              ! add to the general momentum trend  
     520              ua(ji,jj,1) = ua(ji,jj,1) + zhpi(ji,jj,1) + zuap * umask(ji,jj,1)  
     521              va(ji,jj,1) = va(ji,jj,1) + zhpj(ji,jj,1) + zvap * vmask(ji,jj,1)  
     522           END DO  
     523        END DO  
     524  
     525        ! interior value (2=<jk=<jpkm1)  
     526        DO jk = 2, jpkm1  
     527           DO jj = 2, jpjm1  
     528              DO ji = fs_2, fs_jpim1   ! vector opt.  
     529                 ! hydrostatic pressure gradient along s-surfaces  
     530                 zhpi(ji,jj,jk) = zhpi(ji,jj,jk-1) + zcoef0 * r1_e1u(ji,jj)   &  
     531                    &           * (  e3w_n(ji+1,jj,jk) * ( 2*znad )   &  
     532                    &              - e3w_n(ji  ,jj,jk) * ( 2*znad )  )  
     533                 zhpj(ji,jj,jk) = zhpj(ji,jj,jk-1) + zcoef0 * r1_e2v(ji,jj)   &  
     534                    &           * (  e3w_n(ji,jj+1,jk) * ( 2*znad )   &  
     535                    &              - e3w_n(ji,jj  ,jk) * ( 2*znad )  )  
     536                 ! s-coordinate pressure gradient correction  
     537                 zuap = -zcoef0 * ( 2._wp * znad )   &  
     538                    &           * ( gde3w_n(ji+1,jj  ,jk) - gde3w_n(ji,jj,jk) ) * r1_e1u(ji,jj)  
     539                 zvap = -zcoef0 * ( 2._wp * znad )   &  
     540                    &           * ( gde3w_n(ji  ,jj+1,jk) - gde3w_n(ji,jj,jk) ) * r1_e2v(ji,jj)  
     541  
     542                 IF( ln_wd_il ) THEN  
     543                   zhpi(ji,jj,jk) = zhpi(ji,jj,jk) * zcpx(ji,jj)  
     544                   zhpj(ji,jj,jk) = zhpj(ji,jj,jk) * zcpy(ji,jj)   
     545                   zuap = zuap * zcpx(ji,jj)  
     546                   zvap = zvap * zcpy(ji,jj)  
     547                 ENDIF  
     548  
     549                 ! add to the general momentum trend  
     550                 ua(ji,jj,jk) = ua(ji,jj,jk) + zhpi(ji,jj,jk) + zuap * umask(ji,jj,1)  
     551                 va(ji,jj,jk) = va(ji,jj,jk) + zhpj(ji,jj,jk) + zvap * vmask(ji,jj,1)  
     552              END DO  
     553           END DO  
     554        END DO  
     555  
     556    ELSE  
     557 
    495558      ! Surface value 
    496559      DO jj = 2, jpjm1 
     
    550613         END DO 
    551614      END DO 
     615 
     616    ENDIF !ln_2d 
    552617      ! 
    553618      IF( ln_wd_il )  DEALLOCATE( zcpx , zcpy ) 
  • NEMO/branches/UKMO/NEMO_4.0_surge/src/OCE/DYN/dynkeg.F90

    r11081 r11180  
    191191         DO jj = 2, jpjm1 
    192192            DO ji = fs_2, fs_jpim1   ! vector opt. 
    193                ua(ji,jj,jk) = ua(ji,jj,jk) - ( zhke(ji+1,jj  ,jk) - zhke(ji,jj,jk) ) / e1u(ji,jj) 
    194                va(ji,jj,jk) = va(ji,jj,jk) - ( zhke(ji  ,jj+1,jk) - zhke(ji,jj,jk) ) / e2v(ji,jj) 
     193               ua(ji,jj,jk) = ua(ji,jj,jk) - ( zhke(ji+1,jj  ,jk) - zhke(ji,jj,jk) ) / e1u(ji,jj)* umask(ji,jj,jk) 
     194               va(ji,jj,jk) = va(ji,jj,jk) - ( zhke(ji  ,jj+1,jk) - zhke(ji,jj,jk) ) / e2v(ji,jj)* vmask(ji,jj,jk) 
    195195            END DO  
    196196         END DO 
  • NEMO/branches/UKMO/NEMO_4.0_surge/src/OCE/DYN/dynspg.F90

    r10888 r11180  
    150150            DO jj = 2, jpjm1 
    151151               DO ji = fs_2, fs_jpim1   ! vector opt. 
    152                   ua(ji,jj,jk) = ua(ji,jj,jk) + spgu(ji,jj) 
    153                   va(ji,jj,jk) = va(ji,jj,jk) + spgv(ji,jj) 
     152                  ua(ji,jj,jk) = ua(ji,jj,jk) + spgu(ji,jj)* umask(ji,jj,jk) 
     153                  va(ji,jj,jk) = va(ji,jj,jk) + spgv(ji,jj)* vmask(ji,jj,jk) 
    154154               END DO 
    155155            END DO 
  • NEMO/branches/UKMO/NEMO_4.0_surge/src/OCE/SBC/sbc_oce.F90

    r10888 r11180  
    112112   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   vtau   , vtau_b   !: sea surface j-stress (ocean referential)     [N/m2] 
    113113   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   taum              !: module of sea surface stress (at T-point)    [N/m2]  
     114   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   uwnd              !: wind speed i-component                       [m/s]  
     115   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   vwnd              !: wind speed j-component                       [m/s]  
    114116   !! wndm is used onmpute surface gases exchanges in ice-free ocean or leads 
    115117   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   wndm              !: wind speed module at T-point (=|U10m-Uoce|)  [m/s] 
     
    167169      ! 
    168170      ALLOCATE( utau(jpi,jpj) , utau_b(jpi,jpj) , taum(jpi,jpj) ,     & 
    169          &      vtau(jpi,jpj) , vtau_b(jpi,jpj) , wndm(jpi,jpj) , STAT=ierr(1) )  
    170          ! 
     171         &      vtau(jpi,jpj) , vtau_b(jpi,jpj) , wndm(jpi,jpj) ,     &  
     172         &      uwnd(jpi,jpj) , vwnd(jpi,jpj)   , STAT=ierr(1) )   
     173          ! 
    171174      ALLOCATE( qns_tot(jpi,jpj) , qns  (jpi,jpj) , qns_b(jpi,jpj),        & 
    172175         &      qsr_tot(jpi,jpj) , qsr  (jpi,jpj) ,                        & 
  • NEMO/branches/UKMO/NEMO_4.0_surge/src/OCE/SBC/tide.h90

    r10888 r11180  
    11   !!---------------------------------------------------------------------- 
    22   !! History :  3.2  !  2007  (O. Le Galloudec)  Original code 
     3   !!                 !  2019  (C O'Neill) Added components 20-31  
    34   !!---------------------------------------------------------------------- 
    45 
     
    2829   Wave(18) = tide(  'L2'     , 0.006694 ,    2   ,  2 , -1 ,  2 , -1 ,  0  , +180  ,  2   , -2   ,  0   ,  0   , 0 ,  215    ) 
    2930   Wave(19) = tide(  'T2'     , 0.006614 ,    2   ,  2 ,  0 , -1 ,  0 ,  1  ,    0  ,  0   ,  0   ,  0   ,  0   , 0 ,    0    ) 
     31   !  
     32   Wave(20) = tide(  'MNS2'   , 0.000000 ,    2   ,  2 , -5 ,  4 ,  1 ,  0  ,    0  ,  4   , -4   ,  0   ,  0   , 0 ,    6    )  
     33   Wave(21) = tide(  'Lam2'   , 0.001760 ,    2   ,  2 , -1 ,  0 ,  1 ,  0  , +180  ,  2   , -2   ,  0   ,  0   , 0 ,   78    )  
     34   Wave(22) = tide(  'MSN2'   , 0.000000 ,    2   ,  2 ,  1 ,  0 ,  1 ,  0  ,    0  ,  2   , -2   ,  0   ,  2   , 0 ,    6    )  
     35   Wave(23) = tide(  '2SM2'   , 0.000000 ,    2   ,  2 ,  2 , -2 ,  0 ,  0  ,    0  , -2   ,  2   ,  0   ,  0   , 0 ,   16    )  
     36   Wave(24) = tide(  'MO3'    , 0.000000 ,    3   ,  3 , -4 ,  1 ,  0 ,  0  ,  +90  ,  2   , -2   ,  0   ,  0   , 0 ,   13    )  
     37   Wave(25) = tide(  'MK3'    , 0.000000 ,    3   ,  3 , -2 ,  3 ,  0 ,  0  ,  -90  ,  2   , -2   , -1   ,  0   , 0 ,   10    )  
     38   Wave(26) = tide(  'MN4'    , 0.000000 ,    4   ,  4 , -5 ,  4 ,  1 ,  0  ,    0  ,  4   , -4   ,  0   ,  0   , 0 ,    1    )  
     39   Wave(27) = tide(  'MS4'    , 0.000000 ,    4   ,  4 , -2 ,  2 ,  0 ,  0  ,    0  ,  2   , -2   ,  0   ,  0   , 0 ,    2    )  
     40   Wave(28) = tide(  'M6'     , 0.000000 ,    6   ,  6 , -6 ,  6 ,  0 ,  0  ,    0  ,  6   , -6   ,  0   ,  0   , 0 ,    4    )  
     41   Wave(29) = tide(  '2MS6'   , 0.000000 ,    6   ,  6 , -4 ,  4 ,  0 ,  0  ,    0  ,  4   , -4   ,  0   ,  0   , 0 ,    6    )  
     42   Wave(30) = tide(  '2MK6'   , 0.000000 ,    6   ,  6 , -4 ,  6 ,  0 ,  0  ,    0  ,  4   , -4   ,  0   , -2   , 0 ,    5    )  
     43   Wave(31) = tide(  '3M2S2'  , 0.000000 ,    2   , 2  , -6 ,  6 ,  0 ,  0  ,    0  ,  6   , -6   ,  0   ,  0   , 0 ,   12    )  
  • NEMO/branches/UKMO/NEMO_4.0_surge/src/OCE/SBC/tide_mod.F90

    r10888 r11180  
    1616   PUBLIC   tide_init_Wave   ! called by tideini and diaharm modules 
    1717 
    18    INTEGER, PUBLIC, PARAMETER ::   jpmax_harmo = 19   !: maximum number of harmonic 
     18   INTEGER, PUBLIC, PARAMETER ::   jpmax_harmo = 31   !: maximum number of harmonic 
    1919 
    2020   TYPE, PUBLIC ::    tide 
    21       CHARACTER(LEN=4) ::   cname_tide 
     21      CHARACTER(LEN=5) ::   cname_tide 
    2222      REAL(wp)         ::   equitide 
    2323      INTEGER          ::   nutide 
  • NEMO/branches/UKMO/NEMO_4.0_surge/src/OCE/SBC/tideini.F90

    r10888 r11180  
    5050      !!----------------------------------------------------------------------       
    5151      INTEGER  :: ji, jk 
    52       CHARACTER(LEN=4), DIMENSION(jpmax_harmo) :: clname 
     52      CHARACTER(LEN=5), DIMENSION(jpmax_harmo) :: clname 
    5353      INTEGER  ::   ios                 ! Local integer output status for namelist read 
    5454      ! 
  • NEMO/branches/UKMO/NEMO_4.0_surge/src/OCE/ZDF/zdfphy.F90

    r11081 r11180  
    8484      !!---------------------------------------------------------------------- 
    8585      ! 
     86      IF( .NOT. ln_2d ) THEN 
     87 
    8688      IF(lwp) THEN 
    8789         WRITE(numout,*) 
     
    209211      IF( ln_zdfswm )   CALL zdf_swm_init       ! surface  wave-driven mixing 
    210212 
     213      ENDIF  ! .NOT. ln_2d 
     214 
    211215      !                          !== top/bottom friction  ==! 
    212216      CALL zdf_drg_init 
     
    251255         ENDIF 
    252256      ENDIF 
     257       
     258      IF( .NOT. ln_2d ) THEN ! 2D case only uses bottom friction from this routine      
    253259      ! 
    254260      !                       !==  Kz from chosen turbulent closure  ==!   (avm_k, avt_k) 
     
    321327         ! NB. OSMOSIS restart (osm_rst) will be called in step.F90 after wn has been updated 
    322328      ENDIF 
     329       
     330      ENDIF  ! .NOT. ln_2d 
    323331      ! 
    324332      IF( ln_timing )   CALL timing_stop('zdf_phy') 
  • NEMO/branches/UKMO/NEMO_4.0_surge/src/OCE/nemogcm.F90

    r10888 r11180  
    7575   USE sbc_oce , ONLY : lk_oasis 
    7676   USE wet_dry        ! Wetting and drying setting   (wad_init routine) 
     77   USE dom_oce , ONLY : ln_2d  
    7778#if defined key_top 
    7879   USE trcini         ! passive tracer initialisation 
     
    446447                                      
    447448      !                                         ! Lateral physics 
    448                            CALL ldf_tra_init      ! Lateral ocean tracer physics 
     449      IF (.NOT. ln_2d)     CALL ldf_tra_init      ! Lateral ocean tracer physics 
    449450                           CALL ldf_eiv_init      ! eddy induced velocity param. 
    450451                           CALL ldf_dyn_init      ! Lateral ocean momentum physics 
    451452 
    452453      !                                      ! Active tracers 
    453       IF( ln_traqsr    )   CALL tra_qsr_init      ! penetrative solar radiation qsr 
    454                            CALL tra_bbc_init      ! bottom heat flux 
    455                            CALL tra_bbl_init      ! advective (and/or diffusive) bottom boundary layer scheme 
    456                            CALL tra_dmp_init      ! internal tracer damping 
    457                            CALL tra_adv_init      ! horizontal & vertical advection 
    458                            CALL tra_ldf_init      ! lateral mixing 
     454      IF (.NOT. ln_2d) THEN  ! No tracers in 2D mode  
     455          IF( ln_traqsr    )   CALL tra_qsr_init      ! penetrative solar radiation qsr 
     456                               CALL tra_bbc_init      ! bottom heat flux 
     457                               CALL tra_bbl_init      ! advective (and/or diffusive) bottom boundary layer scheme 
     458                               CALL tra_dmp_init      ! internal tracer damping 
     459                               CALL tra_adv_init      ! horizontal & vertical advection 
     460                               CALL tra_ldf_init      ! lateral mixing 
     461      ENDIF 
    459462 
    460463      !                                      ! Dynamics 
  • NEMO/branches/UKMO/NEMO_4.0_surge/src/OCE/step.F90

    r10888 r11180  
    125125      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    126126      !  THERMODYNAMICS 
     127      IF ( .NOT. ln_2d ) THEN  
    127128                         CALL eos_rab( tsb, rab_b )       ! before local thermal/haline expension ratio at T-points 
    128129                         CALL eos_rab( tsn, rab_n )       ! now    local thermal/haline expension ratio at T-points 
    129130                         CALL bn2    ( tsb, rab_b, rn2b ) ! before Brunt-Vaisala frequency 
    130131                         CALL bn2    ( tsn, rab_n, rn2  ) ! now    Brunt-Vaisala frequency 
    131  
     132      ENDIF 
     133       
    132134      !  VERTICAL PHYSICS 
    133135                         CALL zdf_phy( kstp )         ! vertical physics update (top/bot drag, avt, avs, avm + MLD) 
     
    163165                            CALL wzv           ( kstp )  ! now cross-level velocity  
    164166      IF( ln_zad_Aimp )     CALL wAimp         ( kstp )  ! Adaptive-implicit vertical advection partitioning 
    165                             CALL eos    ( tsn, rhd, rhop, gdept_n(:,:,:) )  ! now in situ density for hpg computation 
     167      IF( .NOT. ln_2d )     CALL eos    ( tsn, rhd, rhop, gdept_n(:,:,:) )  ! now in situ density for hpg computation 
    166168                             
    167169!!jc: fs simplification 
     
    202204      ENDIF 
    203205       
    204                          CALL dyn_zdf       ( kstp )  ! vertical diffusion 
     206      IF( .NOT. ln_2d )     CALL dyn_zdf       ( kstp )  ! vertical diffusion 
    205207 
    206208      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
     
    232234      ! Active tracers                               
    233235      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
     236      IF (.NOT. ln_2d) THEN  ! No tracers in 2D mode 
    234237                         tsa(:,:,:,:) = 0._wp         ! set tracer trends to zero 
    235238 
     
    257260                         CALL tra_zdf       ( kstp )  ! vertical mixing and after tracer fields 
    258261      IF( ln_zdfnpc  )   CALL tra_npc       ( kstp )  ! update after fields by non-penetrative convection 
    259  
     262      ENDIF ! not ln_2d  
     263       
    260264      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
    261265      ! Set boundary conditions and Swap 
     
    275279!!  
    276280!!jc2: dynnxt must be the latest call. e3t_b are indeed updated in that routine 
    277                          CALL tra_nxt       ( kstp )  ! finalize (bcs) tracer fields at next time step and swap 
     281      IF (.NOT. ln_2d)   CALL tra_nxt       ( kstp )  ! finalize (bcs) tracer fields at next time step and swap 
    278282                         CALL dyn_nxt       ( kstp )  ! finalize (bcs) velocities at next time step and swap (always called after tra_nxt) 
    279283                         CALL ssh_swp       ( kstp )  ! swap of sea surface height 
Note: See TracChangeset for help on using the changeset viewer.