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 11949 for NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/C1D – NEMO

Ignore:
Timestamp:
2019-11-22T15:29:17+01:00 (5 years ago)
Author:
acc
Message:

Merge in changes from 2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps. This just creates a fresh copy of this branch to use as the merge base. See ticket #2341

Location:
NEMO/branches/2019/dev_r11943_MERGE_2019/src
Files:
5 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2019/dev_r11943_MERGE_2019/src

    • Property svn:mergeinfo deleted
  • NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/C1D/dtauvd.F90

    r11536 r11949  
    117117 
    118118 
    119    SUBROUTINE dta_uvd( kt, puvd ) 
     119   SUBROUTINE dta_uvd( kt, Kmm, puvd ) 
    120120      !!---------------------------------------------------------------------- 
    121121      !!                   ***  ROUTINE dta_uvd  *** 
     
    133133      !!---------------------------------------------------------------------- 
    134134      INTEGER                           , INTENT(in   ) ::   kt     ! ocean time-step 
     135      INTEGER                           , INTENT(in   ) ::   Kmm    ! time level index 
    135136      REAL(wp), DIMENSION(jpi,jpj,jpk,2), INTENT(  out) ::   puvd   ! U & V current data 
    136137      ! 
     
    160161            DO ji = 1, jpi                ! determines the interpolated U & V current profiles at each (i,j) point 
    161162               DO jk = 1, jpk 
    162                   zl = gdept_n(ji,jj,jk) 
     163                  zl = gdept(ji,jj,jk,Kmm) 
    163164                  IF    ( zl < gdept_1d(1  ) ) THEN          ! extrapolate above the first level of data 
    164165                     zup(jk) =  puvd(ji,jj,1    ,1) 
  • NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/C1D/dyncor_c1d.F90

    r10068 r11949  
    5656 
    5757 
    58    SUBROUTINE dyn_cor_c1d( kt ) 
     58   SUBROUTINE dyn_cor_c1d( kt, Kmm, puu, pvv, Krhs ) 
    5959      !!---------------------------------------------------------------------- 
    6060      !!                   ***  ROUTINE dyn_cor_c1d  *** 
     
    6363      !!               the general trend of the momentum equation in 1D case. 
    6464      !!---------------------------------------------------------------------- 
    65       INTEGER, INTENT( in ) ::   kt   ! ocean time-step index 
     65      INTEGER                             , INTENT(in   ) ::   kt        ! ocean time-step index 
     66      INTEGER                             , INTENT(in   ) ::   Kmm, Krhs ! ocean time level indices 
     67      REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT(inout) ::   puu, pvv  ! ocean velocities and RHS of momentum equation 
    6668      !! 
    6769      INTEGER ::   ji, jj, jk   ! dummy loop indices 
     
    7880            DO jj = 2, jpjm1 
    7981               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                  puu(ji,jj,jk,Krhs) = puu(ji,jj,jk,Krhs) + ff_t(ji,jj) * (pvv(ji,jj,jk,Kmm) + vsd(ji,jj,jk)) 
     83                  pvv(ji,jj,jk,Krhs) = pvv(ji,jj,jk,Krhs) - ff_t(ji,jj) * (puu(ji,jj,jk,Kmm) + usd(ji,jj,jk)) 
    8284               END DO 
    8385            END DO 
     
    8789            DO jj = 2, jpjm1 
    8890               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                  puu(ji,jj,jk,Krhs) = puu(ji,jj,jk,Krhs) + ff_t(ji,jj) * pvv(ji,jj,jk,Kmm) 
     92                  pvv(ji,jj,jk,Krhs) = pvv(ji,jj,jk,Krhs) - ff_t(ji,jj) * puu(ji,jj,jk,Kmm) 
    9193               END DO 
    9294            END DO 
     
    9597       
    9698      ! 
    97       IF(ln_ctl)   CALL prt_ctl( tab3d_1=ua, clinfo1=' cor  - Ua: ', mask1=umask,  & 
    98          &                       tab3d_2=va, clinfo2=' Va: '       , mask2=vmask ) 
     99      IF(ln_ctl)   CALL prt_ctl( tab3d_1=puu(:,:,:,Krhs), clinfo1=' cor  - Ua: ', mask1=umask,  & 
     100         &                       tab3d_2=pvv(:,:,:,Krhs), clinfo2=' Va: '       , mask2=vmask ) 
    99101      ! 
    100102   END SUBROUTINE dyn_cor_c1d 
  • NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/C1D/dyndmp.F90

    r11536 r11949  
    130130 
    131131 
    132    SUBROUTINE dyn_dmp( kt ) 
     132   SUBROUTINE dyn_dmp( kt, Kbb, Kmm, puu, pvv, Krhs ) 
    133133      !!---------------------------------------------------------------------- 
    134134      !!                   ***  ROUTINE dyn_dmp  *** 
     
    140140      !! ** Method  :   Compute Newtonian damping towards u_dta and v_dta  
    141141      !!      and add to the general momentum trends: 
    142       !!                     ua = ua + resto_uv * (u_dta - ub) 
    143       !!                     va = va + resto_uv * (v_dta - vb) 
     142      !!                     puu(Krhs) = puu(Krhs) + resto_uv * (u_dta - puu(Kbb)) 
     143      !!                     pvv(Krhs) = pvv(Krhs) + resto_uv * (v_dta - pvv(Kbb)) 
    144144      !!      The trend is computed either throughout the water column 
    145145      !!      (nn_zdmp=0), where the vertical mixing is weak (nn_zdmp=1) or 
    146146      !!      below the well mixed layer (nn_zdmp=2) 
    147147      !! 
    148       !! ** Action  : - (ua,va)   momentum trends updated with the damping trend 
    149       !!---------------------------------------------------------------------- 
    150       INTEGER, INTENT(in) ::   kt   ! ocean time-step index 
     148      !! ** Action  : - (puu(:,:,:,Krhs),pvv(:,:,:,Krhs))   momentum trends updated with the damping trend 
     149      !!---------------------------------------------------------------------- 
     150      INTEGER                             , INTENT(in   ) ::   kt             ! ocean time-step index 
     151      INTEGER                             , INTENT(in   ) ::   Kbb, Kmm, Krhs ! ocean time level indices 
     152      REAL(wp), DIMENSION(jpi,jpj,jpk,jpt), INTENT(inout) ::   puu, pvv       ! ocean velocities and RHS of momentum equation 
    151153      !! 
    152154      INTEGER  ::   ji, jj, jk   ! dummy loop indices 
     
    159161      ! 
    160162      !                           !==   read and interpolate U & V current data at kt   ==! 
    161       CALL dta_uvd( kt, zuv_dta ) !!! NOTE: This subroutine must be altered for use outside 
     163      CALL dta_uvd( kt, Kmm, zuv_dta ) !!! NOTE: This subroutine must be altered for use outside 
    162164                                  !!!       the C1D context (use of U,V grid variables) 
    163165      ! 
     
    168170            DO jj = 2, jpjm1 
    169171               DO ji = fs_2, fs_jpim1   ! vector opt. 
    170                   zua = resto_uv(ji,jj,jk) * ( zuv_dta(ji,jj,jk,1) - ub(ji,jj,jk) ) 
    171                   zva = resto_uv(ji,jj,jk) * ( zuv_dta(ji,jj,jk,2) - vb(ji,jj,jk) ) 
    172                   ua(ji,jj,jk) = ua(ji,jj,jk) + zua 
    173                   va(ji,jj,jk) = va(ji,jj,jk) + zva 
     172                  zua = resto_uv(ji,jj,jk) * ( zuv_dta(ji,jj,jk,1) - puu(ji,jj,jk,Kbb) ) 
     173                  zva = resto_uv(ji,jj,jk) * ( zuv_dta(ji,jj,jk,2) - pvv(ji,jj,jk,Kbb) ) 
     174                  puu(ji,jj,jk,Krhs) = puu(ji,jj,jk,Krhs) + zua 
     175                  pvv(ji,jj,jk,Krhs) = pvv(ji,jj,jk,Krhs) + zva 
    174176                  utrdmp(ji,jj,jk) = zua           ! save the trends 
    175177                  vtrdmp(ji,jj,jk) = zva       
     
    183185               DO ji = fs_2, fs_jpim1   ! vector opt. 
    184186                  IF( avt(ji,jj,jk) <= avt_c ) THEN 
    185                      zua = resto_uv(ji,jj,jk) * ( zuv_dta(ji,jj,jk,1) - ub(ji,jj,jk) ) 
    186                      zva = resto_uv(ji,jj,jk) * ( zuv_dta(ji,jj,jk,2) - vb(ji,jj,jk) ) 
     187                     zua = resto_uv(ji,jj,jk) * ( zuv_dta(ji,jj,jk,1) - puu(ji,jj,jk,Kbb) ) 
     188                     zva = resto_uv(ji,jj,jk) * ( zuv_dta(ji,jj,jk,2) - pvv(ji,jj,jk,Kbb) ) 
    187189                  ELSE 
    188190                     zua = 0._wp 
    189191                     zva = 0._wp   
    190192                  ENDIF 
    191                   ua(ji,jj,jk) = ua(ji,jj,jk) + zua 
    192                   va(ji,jj,jk) = va(ji,jj,jk) + zva 
     193                  puu(ji,jj,jk,Krhs) = puu(ji,jj,jk,Krhs) + zua 
     194                  pvv(ji,jj,jk,Krhs) = pvv(ji,jj,jk,Krhs) + zva 
    193195                  utrdmp(ji,jj,jk) = zua           ! save the trends 
    194196                  vtrdmp(ji,jj,jk) = zva 
     
    201203            DO jj = 2, jpjm1 
    202204               DO ji = fs_2, fs_jpim1   ! vector opt. 
    203                   IF( gdept_n(ji,jj,jk) >= hmlp (ji,jj) ) THEN 
    204                      zua = resto_uv(ji,jj,jk) * ( zuv_dta(ji,jj,jk,1) - ub(ji,jj,jk) ) 
    205                      zva = resto_uv(ji,jj,jk) * ( zuv_dta(ji,jj,jk,2) - vb(ji,jj,jk) ) 
     205                  IF( gdept(ji,jj,jk,Kmm) >= hmlp (ji,jj) ) THEN 
     206                     zua = resto_uv(ji,jj,jk) * ( zuv_dta(ji,jj,jk,1) - puu(ji,jj,jk,Kbb) ) 
     207                     zva = resto_uv(ji,jj,jk) * ( zuv_dta(ji,jj,jk,2) - pvv(ji,jj,jk,Kbb) ) 
    206208                  ELSE 
    207209                     zua = 0._wp 
    208210                     zva = 0._wp   
    209211                  ENDIF 
    210                   ua(ji,jj,jk) = ua(ji,jj,jk) + zua 
    211                   va(ji,jj,jk) = va(ji,jj,jk) + zva 
     212                  puu(ji,jj,jk,Krhs) = puu(ji,jj,jk,Krhs) + zua 
     213                  pvv(ji,jj,jk,Krhs) = pvv(ji,jj,jk,Krhs) + zva 
    212214                  utrdmp(ji,jj,jk) = zua           ! save the trends 
    213215                  vtrdmp(ji,jj,jk) = zva 
     
    219221      ! 
    220222      !                           ! Control print 
    221       IF( ln_ctl   )   CALL prt_ctl( tab3d_1=ua(:,:,:), clinfo1=' dmp  - Ua: ', mask1=umask,   & 
    222          &                           tab3d_2=va(:,:,:), clinfo2=       ' Va: ', mask2=vmask, clinfo3='dyn' ) 
     223      IF( ln_ctl   )   CALL prt_ctl( tab3d_1=puu(:,:,:,Krhs), clinfo1=' dmp  - Ua: ', mask1=umask,   & 
     224         &                           tab3d_2=pvv(:,:,:,Krhs), clinfo2=       ' Va: ', mask2=vmask, clinfo3='dyn' ) 
    223225      ! 
    224226      ! 
  • NEMO/branches/2019/dev_r11943_MERGE_2019/src/OCE/C1D/step_c1d.F90

    r10068 r11949  
    66   !! History :   2.0  !  2004-04  (C. Ethe)  adapted from step.F90 for C1D 
    77   !!             3.0  !  2008-04  (G. Madec)  redo the adaptation to include SBC 
     8   !!             4.1  !  2019-08  (A. Coward, D. Storkey) rewrite in preparation for new timestepping scheme 
    89   !!---------------------------------------------------------------------- 
    910#if defined key_c1d 
     
    1415   !!---------------------------------------------------------------------- 
    1516   USE step_oce        ! time stepping definition modules  
     17   USE step, ONLY : Nbb, Nnn, Naa, Nrhs ! time level indices 
    1618#if defined key_top 
    1719   USE trcstp          ! passive tracer time-stepping      (trc_stp routine) 
    1820#endif 
    1921   USE dyncor_c1d      ! Coriolis term (c1d case)         (dyn_cor_1d     ) 
    20    USE dynnxt          ! time-stepping                    (dyn_nxt routine) 
     22   USE dynatf          ! time filtering                   (dyn_atf routine) 
    2123   USE dyndmp          ! U & V momentum damping           (dyn_dmp routine) 
    2224   USE restart         ! restart  
     
    6567      ! Update data, open boundaries, surface boundary condition (including sea-ice) 
    6668      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    67                          CALL sbc    ( kstp )         ! Sea Boundary Condition (including sea-ice) 
     69                         CALL sbc    ( kstp, Nbb, Nnn )  ! Sea Boundary Condition (including sea-ice) 
    6870 
    6971      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
    70       ! Ocean physics update                (ua, va, ta, sa used as workspace) 
     72      ! Ocean physics update         
    7173      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    72                          CALL eos_rab( tsb, rab_b )   ! before local thermal/haline expension ratio at T-points 
    73                          CALL eos_rab( tsn, rab_n )   ! now    local thermal/haline expension ratio at T-points 
    74                          CALL bn2( tsb, rab_b, rn2b ) ! before Brunt-Vaisala frequency 
    75                          CALL bn2( tsn, rab_n, rn2 ) ! now    Brunt-Vaisala frequency 
     74                         CALL eos_rab( ts(:,:,:,:,Nbb), rab_b, Nnn )  ! before local thermal/haline expension ratio at T-points 
     75                         CALL eos_rab( ts(:,:,:,:,Nnn), rab_n, Nnn )  ! now    local thermal/haline expension ratio at T-points 
     76                         CALL bn2( ts(:,:,:,:,Nbb), rab_b, rn2b, Nnn ) ! before Brunt-Vaisala frequency 
     77                         CALL bn2( ts(:,:,:,:,Nnn), rab_n, rn2 , Nnn ) ! now    Brunt-Vaisala frequency 
    7678       
    7779      !  VERTICAL PHYSICS 
    78                          CALL zdf_phy( kstp )         ! vertical physics update (bfr, avt, avs, avm + MLD) 
     80                         CALL zdf_phy( kstp, Nbb, Nnn, Nrhs  )    ! vertical physics update (bfr, avt, avs, avm + MLD) 
    7981 
    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      IF(.NOT.ln_linssh )   CALL ssh_nxt       ( kstp, Nbb, Nnn, ssh, Naa )  ! after ssh (includes call to div_hor) 
     83      IF(.NOT.ln_linssh )   CALL dom_vvl_sf_nxt( kstp, Nbb, Nnn,      Naa )  ! after vertical scale factors  
    8284 
    83       IF(.NOT.ln_linssh )   CALL wzv           ( kstp )  ! now cross-level velocity  
     85      IF(.NOT.ln_linssh )   CALL wzv           ( kstp, Nbb, Nnn, ww,  Naa )  ! now cross-level velocity  
    8486      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
    85       ! diagnostics and outputs             (ua, va, ta, sa used as workspace) 
     87      ! diagnostics and outputs        
    8688      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    87                          CALL dia_wri( kstp )       ! ocean model: outputs 
    88       IF( lk_diahth  )   CALL dia_hth( kstp )       ! Thermocline depth (20°C) 
     89                         CALL dia_wri( kstp, Nnn )  ! ocean model: outputs 
     90      IF( lk_diahth  )   CALL dia_hth( kstp, Nnn )  ! Thermocline depth (20°C) 
    8991 
    9092 
     
    9395      ! Passive Tracer Model 
    9496      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    95                         CALL trc_stp( kstp )       ! time-stepping 
     97                        CALL trc_stp( kstp, Nbb, Nnn, Nrhs, Naa  )   ! time-stepping 
    9698#endif 
    9799 
    98100      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
    99       ! Active tracers                              (ua, va used as workspace) 
     101      ! Active tracers                              (uu(:,:,:,Nrhs), vv(:,:,:,Nrhs) used as workspace) 
    100102      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    101                         tsa(:,:,:,:) = 0._wp       ! set tracer trends to zero 
     103                        ts(:,:,:,:,Nrhs) = 0._wp       ! set tracer trends to zero 
    102104 
    103                         CALL tra_sbc( kstp )       ! surface boundary condition 
    104       IF( ln_traqsr )   CALL tra_qsr( kstp )       ! penetrative solar radiation qsr 
    105       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( ln_zdfosm  )  CALL tra_osm( kstp )       ! OSMOSIS non-local tracer fluxes 
    108                         CALL tra_zdf( kstp )       ! vertical mixing 
    109                         CALL eos( tsn, rhd, rhop, gdept_0(:,:,:) )   ! now potential density for zdfmxl 
    110       IF( ln_zdfnpc )   CALL tra_npc( kstp )       ! applied non penetrative convective adjustment on (t,s) 
    111                         CALL tra_nxt( kstp )       ! tracer fields at next time step 
     105                        CALL tra_sbc( kstp,      Nnn, ts, Nrhs  )  ! surface boundary condition 
     106      IF( ln_traqsr )   CALL tra_qsr( kstp,      Nnn, ts, Nrhs  )  ! penetrative solar radiation qsr 
     107      IF( ln_tradmp )   CALL tra_dmp( kstp, Nbb, Nnn, ts, Nrhs  )  ! internal damping trends- tracers 
     108      IF(.NOT.ln_linssh)CALL tra_adv( kstp, Nbb, Nnn, ts, Nrhs  )  ! horizontal & vertical advection 
     109      IF( ln_zdfosm  )  CALL tra_osm( kstp, Nnn     , ts, Nrhs  )  ! OSMOSIS non-local tracer fluxes 
     110                        CALL tra_zdf( kstp, Nbb, Nnn, Nrhs, ts, Naa   )         ! vertical mixing 
     111                        CALL eos( ts(:,:,:,:,Nnn), rhd, rhop, gdept_0(:,:,:) )  ! now potential density for zdfmxl 
     112      IF( ln_zdfnpc )   CALL tra_npc( kstp,      Nnn, Nrhs, ts, Naa   )         ! applied non penetrative convective adjustment on (t,s) 
     113                        CALL tra_atf( kstp, Nbb, Nnn, Nrhs,     Naa, ts   )     ! time filtering of "now" tracer fields 
    112114 
    113115 
     
    116118      ! Dynamics                                    (ta, sa used as workspace) 
    117119      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    118                         ua(:,:,:) = 0._wp          ! set dynamics trends to zero 
    119                         va(:,:,:) = 0._wp 
     120                        uu(:,:,:,Nrhs) = 0._wp          ! set dynamics trends to zero 
     121                        vv(:,:,:,Nrhs) = 0._wp 
    120122 
    121       IF( ln_dyndmp )   CALL dyn_dmp    ( kstp )   ! internal damping trends- momentum 
    122                         CALL dyn_cor_c1d( kstp )   ! vorticity term including Coriolis 
    123       IF( ln_zdfosm  )  CALL dyn_osm    ( kstp )   ! OSMOSIS non-local velocity fluxes 
    124                         CALL dyn_zdf    ( kstp )   ! vertical diffusion 
    125                         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 
     123      IF( ln_dyndmp )   CALL dyn_dmp    ( kstp, Nbb, Nnn      , uu, vv, Nrhs )  ! internal damping trends- momentum 
     124                        CALL dyn_cor_c1d( kstp,      Nnn      , uu, vv, Nrhs )  ! vorticity term including Coriolis 
     125      IF( ln_zdfosm  )  CALL dyn_osm    ( kstp,      Nnn      , uu, vv, Nrhs )  ! OSMOSIS non-local velocity fluxes 
     126                        CALL dyn_zdf    ( kstp, Nbb, Nnn, Nrhs, uu, vv, Naa  )  ! vertical diffusion 
     127                        CALL dyn_atf    ( kstp, Nbb, Nnn, Naa , uu, vv, e3t, e3u, e3v )  ! time filtering of "now" fields 
     128      IF(.NOT.ln_linssh)CALL ssh_atf    ( kstp, Nbb, Nnn, Naa , ssh )                    ! time filtering of "now" sea surface height 
     129      ! 
     130      ! Swap time levels 
     131      Nrhs = Nbb 
     132      Nbb = Nnn 
     133      Nnn = Naa 
     134      Naa = Nrhs 
     135      ! 
     136      IF(.NOT.ln_linssh)CALL dom_vvl_sf_update( kstp, Nbb, Nnn, Naa )                    ! update of vertical scale factors 
    129137 
    130138      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
    131139      ! Control and restarts 
    132140      !<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<<< 
    133                              CALL stp_ctl( kstp, indic ) 
    134       IF( kstp == nit000 )   CALL iom_close( numror )      ! close input  ocean restart file 
    135       IF( lrst_oce       )   CALL rst_write( kstp )        ! write output ocean restart file 
     141                             CALL stp_ctl( kstp, Nnn, indic ) 
     142      IF( kstp == nit000 )   CALL iom_close( numror )          ! close input  ocean restart file 
     143      IF( lrst_oce       )   CALL rst_write( kstp, Nbb, Nnn )  ! write output ocean restart file 
    136144      ! 
    137145#if defined key_iomput 
Note: See TracChangeset for help on using the changeset viewer.