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 15249 – NEMO

Changeset 15249


Ignore:
Timestamp:
2021-09-13T11:59:09+02:00 (3 years ago)
Author:
hadcv
Message:

#2721: Fix SETTE debug failures

Location:
NEMO/trunk
Files:
8 edited

Legend:

Unmodified
Added
Removed
  • NEMO/trunk/cfgs/SHARED/field_def_nemo-oce.xml

    r15017 r15249  
    12431243 
    12441244  <!-- 25h diagnostic output --> 
    1245   <field_group id="25h_grid_T" grid_ref="grid_T_3D" operation="instant"> 
     1245  <field_group id="25h_grid_T" grid_ref="grid_T_3D_inner" operation="instant"> 
    12461246    <field id="temper25h"         name="potential temperature 25h mean"    unit="degC" /> 
    12471247    <field id="tempis25h"         name="insitu temperature 25h mean"    unit="degC" /> 
    12481248    <field id="salin25h"          name="salinity 25h mean"                 unit="psu"  /> 
    1249     <field id="ssh25h"            name="sea surface height 25h mean"  grid_ref="grid_T_2D"      unit="m"    /> 
    1250   </field_group> 
    1251  
    1252   <field_group id="25h_grid_U" grid_ref="grid_U_3D" operation="instant" > 
     1249    <field id="ssh25h"            name="sea surface height 25h mean"  grid_ref="grid_T_2D_inner"      unit="m"    /> 
     1250  </field_group> 
     1251 
     1252  <field_group id="25h_grid_U" grid_ref="grid_U_3D_inner" operation="instant" > 
    12531253    <field id="vozocrtx25h"         name="i current 25h mean"    unit="m/s"   /> 
    12541254  </field_group> 
    12551255 
    1256   <field_group id="25h_grid_V" grid_ref="grid_V_3D" operation="instant"> 
     1256  <field_group id="25h_grid_V" grid_ref="grid_V_3D_inner" operation="instant"> 
    12571257    <field id="vomecrty25h"         name="j current 25h mean"    unit="m/s"    /> 
    12581258  </field_group> 
    12591259 
    1260   <field_group id="25h_grid_W" grid_ref="grid_W_3D" operation="instant"> 
     1260  <field_group id="25h_grid_W" grid_ref="grid_W_3D_inner" operation="instant"> 
    12611261    <field id="vovecrtz25h"         name="k current 25h mean"                 unit="m/s"      /> 
    12621262    <field id="avt25h"              name="vertical diffusivity25h mean"       unit="m2/s" /> 
  • NEMO/trunk/cfgs/SHARED/grid_def_nemo.xml

    r15017 r15249  
    8787  </grid> 
    8888  <grid id="grid_W_3D_inner" > 
    89     <domain domain_ref="grid_W" /> 
    90     <axis axis_ref="depthw_inner" /> 
     89    <domain domain_ref="grid_W_inner" /> 
     90    <axis axis_ref="depthw" /> 
    9191  </grid> 
    9292  <!--  --> 
  • NEMO/trunk/src/OCE/DIA/dia25h.F90

    r12489 r15249  
    3232   REAL(wp), SAVE, ALLOCATABLE, DIMENSION(:,:,:) ::   en_25h  , rmxln_25h 
    3333 
     34!! * Substitutions 
     35#  include "do_loop_substitute.h90" 
    3436   !!---------------------------------------------------------------------- 
    3537   !! NEMO/OCE 4.0 , NEMO Consortium (2018) 
     
    5153      INTEGER ::   ios                 ! Local integer output status for namelist read 
    5254      INTEGER ::   ierror              ! Local integer for memory allocation 
     55      INTEGER ::   ji, jj, jk 
    5356      ! 
    5457      NAMELIST/nam_dia25h/ ln_dia25h 
     
    7376      ! ------------------- ! 
    7477      !                                ! ocean arrays 
    75       ALLOCATE( tn_25h (jpi,jpj,jpk), sn_25h (jpi,jpj,jpk), sshn_25h(jpi,jpj)  ,     & 
    76          &      un_25h (jpi,jpj,jpk), vn_25h (jpi,jpj,jpk), wn_25h(jpi,jpj,jpk),     & 
    77          &      avt_25h(jpi,jpj,jpk), avm_25h(jpi,jpj,jpk),                      STAT=ierror ) 
     78      ALLOCATE( tn_25h (A2D(0),jpk), sn_25h (A2D(0),jpk), sshn_25h(A2D(0))  ,     & 
     79         &      un_25h (A2D(0),jpk), vn_25h (A2D(0),jpk), wn_25h(A2D(0),jpk),     & 
     80         &      avt_25h(A2D(0),jpk), avm_25h(A2D(0),jpk),                      STAT=ierror ) 
    7881      IF( ierror > 0 ) THEN 
    7982         CALL ctl_stop( 'dia_25h: unable to allocate ocean arrays' )   ;   RETURN 
    8083      ENDIF 
    8184      IF( ln_zdftke ) THEN             ! TKE physics 
    82          ALLOCATE( en_25h(jpi,jpj,jpk), STAT=ierror ) 
     85         ALLOCATE( en_25h(A2D(0),jpk), STAT=ierror ) 
    8386         IF( ierror > 0 ) THEN 
    8487            CALL ctl_stop( 'dia_25h: unable to allocate en_25h' )   ;   RETURN 
     
    8689      ENDIF 
    8790      IF( ln_zdfgls ) THEN             ! GLS physics 
    88          ALLOCATE( en_25h(jpi,jpj,jpk), rmxln_25h(jpi,jpj,jpk), STAT=ierror ) 
     91         ALLOCATE( en_25h(A2D(0),jpk), rmxln_25h(A2D(0),jpk), STAT=ierror ) 
    8992         IF( ierror > 0 ) THEN 
    9093            CALL ctl_stop( 'dia_25h: unable to allocate en_25h and rmxln_25h' )   ;   RETURN 
     
    9497      ! 2 - Assign Initial Values ! 
    9598      ! ------------------------- ! 
    96       cnt_25h = 1  ! sets the first value of sum at timestep 1 (note - should strictly be at timestep zero so before values used where possible)  
    97       tn_25h  (:,:,:) = ts (:,:,:,jp_tem,Kbb) 
    98       sn_25h  (:,:,:) = ts (:,:,:,jp_sal,Kbb) 
    99       sshn_25h(:,:)   = ssh(:,:,Kbb) 
    100       un_25h  (:,:,:) = uu  (:,:,:,Kbb) 
    101       vn_25h  (:,:,:) = vv  (:,:,:,Kbb) 
    102       avt_25h (:,:,:) = avt (:,:,:) 
    103       avm_25h (:,:,:) = avm (:,:,:) 
     99      cnt_25h = 1  ! sets the first value of sum at timestep 1 (note - should strictly be at timestep zero so before values used where possible) 
     100      DO_3D( 0, 0, 0, 0, 1, jpk ) 
     101         tn_25h (ji,jj,jk) = ts (ji,jj,jk,jp_tem,Kbb) 
     102         sn_25h (ji,jj,jk) = ts (ji,jj,jk,jp_sal,Kbb) 
     103         un_25h (ji,jj,jk) = uu (ji,jj,jk,Kbb) 
     104         vn_25h (ji,jj,jk) = vv (ji,jj,jk,Kbb) 
     105         avt_25h(ji,jj,jk) = avt(ji,jj,jk) 
     106         avm_25h(ji,jj,jk) = avm(ji,jj,jk) 
     107      END_3D 
     108      DO_2D( 0, 0, 0, 0 ) 
     109         sshn_25h(ji,jj) = ssh(ji,jj,Kbb) 
     110      END_2D 
    104111      IF( ln_zdftke ) THEN 
    105          en_25h(:,:,:) = en(:,:,:) 
     112         DO_3D( 0, 0, 0, 0, 1, jpk ) 
     113            en_25h(ji,jj,jk) = en(ji,jj,jk) 
     114         END_3D 
    106115      ENDIF 
    107116      IF( ln_zdfgls ) THEN 
    108          en_25h   (:,:,:) = en    (:,:,:) 
    109          rmxln_25h(:,:,:) = hmxl_n(:,:,:) 
     117         DO_3D( 0, 0, 0, 0, 1, jpk ) 
     118            en_25h   (ji,jj,jk) = en    (ji,jj,jk) 
     119            rmxln_25h(ji,jj,jk) = hmxl_n(ji,jj,jk) 
     120         END_3D 
    110121      ENDIF 
    111122#if defined key_si3 
     
    132143      REAL(wp)                         ::   zsto, zout, zmax, zjulian, zmdi   ! local scalars 
    133144      INTEGER                          ::   i_steps                           ! no of timesteps per hour 
    134       REAL(wp), DIMENSION(jpi,jpj    ) ::   zw2d, un_dm, vn_dm                ! workspace 
    135       REAL(wp), DIMENSION(jpi,jpj,jpk) ::   zw3d                              ! workspace 
    136       REAL(wp), DIMENSION(jpi,jpj,3)   ::   zwtmb                             ! workspace 
     145      REAL(wp), DIMENSION(A2D(0)    ) ::   zw2d, un_dm, vn_dm                ! workspace 
     146      REAL(wp), DIMENSION(A2D(0),jpk) ::   zw3d                              ! workspace 
     147      REAL(wp), DIMENSION(A2D(0),3)    ::   zwtmb                             ! workspace 
    137148      !!---------------------------------------------------------------------- 
    138149 
     
    151162      ! wn_25h could not be initialised in dia_25h_init, so we do it here instead 
    152163      IF( kt == nn_it000 ) THEN 
    153          wn_25h(:,:,:) = ww(:,:,:) 
     164         DO_3D( 0, 0, 0, 0, 1, jpk ) 
     165            wn_25h(ji,jj,jk) = ww(ji,jj,jk) 
     166         END_3D 
    154167      ENDIF 
    155168 
     
    162175         ENDIF 
    163176 
    164          tn_25h  (:,:,:)     = tn_25h  (:,:,:) + ts (:,:,:,jp_tem,Kmm) 
    165          sn_25h  (:,:,:)     = sn_25h  (:,:,:) + ts (:,:,:,jp_sal,Kmm) 
    166          sshn_25h(:,:)       = sshn_25h(:,:)   + ssh(:,:,Kmm) 
    167          un_25h  (:,:,:)     = un_25h  (:,:,:) + uu  (:,:,:,Kmm) 
    168          vn_25h  (:,:,:)     = vn_25h  (:,:,:) + vv  (:,:,:,Kmm) 
    169          wn_25h  (:,:,:)     = wn_25h  (:,:,:) + ww  (:,:,:) 
    170          avt_25h (:,:,:)     = avt_25h (:,:,:) + avt (:,:,:) 
    171          avm_25h (:,:,:)     = avm_25h (:,:,:) + avm (:,:,:) 
     177         DO_3D( 0, 0, 0, 0, 1, jpk ) 
     178            tn_25h  (ji,jj,jk) = tn_25h  (ji,jj,jk) + ts (ji,jj,jk,jp_tem,Kmm) 
     179            sn_25h  (ji,jj,jk) = sn_25h  (ji,jj,jk) + ts (ji,jj,jk,jp_sal,Kmm) 
     180            un_25h  (ji,jj,jk) = un_25h  (ji,jj,jk) + uu (ji,jj,jk,Kmm) 
     181            vn_25h  (ji,jj,jk) = vn_25h  (ji,jj,jk) + vv (ji,jj,jk,Kmm) 
     182            wn_25h  (ji,jj,jk) = wn_25h  (ji,jj,jk) + ww (ji,jj,jk) 
     183            avt_25h (ji,jj,jk) = avt_25h (ji,jj,jk) + avt(ji,jj,jk) 
     184            avm_25h (ji,jj,jk) = avm_25h (ji,jj,jk) + avm(ji,jj,jk) 
     185         END_3D 
     186         DO_2D( 0, 0, 0, 0 ) 
     187            sshn_25h(ji,jj)    = sshn_25h(ji,jj)    + ssh(ji,jj,Kmm) 
     188         END_2D 
    172189         IF( ln_zdftke ) THEN 
    173             en_25h(:,:,:)    = en_25h  (:,:,:) + en(:,:,:) 
     190            DO_3D( 0, 0, 0, 0, 1, jpk ) 
     191               en_25h(ji,jj,jk) = en_25h(ji,jj,jk) + en(ji,jj,jk) 
     192            END_3D 
    174193         ENDIF 
    175194         IF( ln_zdfgls ) THEN 
    176             en_25h   (:,:,:) = en_25h   (:,:,:) + en    (:,:,:) 
    177             rmxln_25h(:,:,:) = rmxln_25h(:,:,:) + hmxl_n(:,:,:) 
     195            DO_3D( 0, 0, 0, 0, 1, jpk ) 
     196               en_25h   (ji,jj,jk) = en_25h   (ji,jj,jk) + en    (ji,jj,jk) 
     197               rmxln_25h(ji,jj,jk) = rmxln_25h(ji,jj,jk) + hmxl_n(ji,jj,jk) 
     198            END_3D 
    178199         ENDIF 
    179200         cnt_25h = cnt_25h + 1 
     
    212233         zmdi=1.e+20 !missing data indicator for masking 
    213234         ! write tracers (instantaneous) 
    214          zw3d(:,:,:) = tn_25h(:,:,:)*tmask(:,:,:) + zmdi*(1.0-tmask(:,:,:)) 
     235         DO_3D( 0, 0, 0, 0, 1, jpk ) 
     236            zw3d(ji,jj,jk) = tn_25h(ji,jj,jk)*tmask(ji,jj,jk) + zmdi*(1.0-tmask(ji,jj,jk)) 
     237         END_3D 
    215238         CALL iom_put("temper25h", zw3d)   ! potential temperature 
    216          zw3d(:,:,:) = sn_25h(:,:,:)*tmask(:,:,:) + zmdi*(1.0-tmask(:,:,:)) 
     239         DO_3D( 0, 0, 0, 0, 1, jpk ) 
     240            zw3d(ji,jj,jk) = sn_25h(ji,jj,jk)*tmask(ji,jj,jk) + zmdi*(1.0-tmask(ji,jj,jk)) 
     241         END_3D 
    217242         CALL iom_put( "salin25h", zw3d  )   ! salinity 
    218          zw2d(:,:) = sshn_25h(:,:)*tmask(:,:,1) + zmdi*(1.0-tmask(:,:,1)) 
     243         DO_2D( 0, 0, 0, 0 ) 
     244            zw2d(ji,jj) = sshn_25h(ji,jj)*tmask(ji,jj,1) + zmdi*(1.0-tmask(ji,jj,1)) 
     245         END_2D 
    219246         IF( ll_wd ) THEN 
    220247            CALL iom_put( "ssh25h", zw2d+ssh_ref )   ! sea surface  
     
    223250         ENDIF 
    224251         ! Write velocities (instantaneous) 
    225          zw3d(:,:,:) = un_25h(:,:,:)*umask(:,:,:) + zmdi*(1.0-umask(:,:,:)) 
     252         DO_3D( 0, 0, 0, 0, 1, jpk ) 
     253            zw3d(ji,jj,jk) = un_25h(ji,jj,jk)*umask(ji,jj,jk) + zmdi*(1.0-umask(ji,jj,jk)) 
     254         END_3D 
    226255         CALL iom_put("vozocrtx25h", zw3d)    ! i-current 
    227          zw3d(:,:,:) = vn_25h(:,:,:)*vmask(:,:,:) + zmdi*(1.0-vmask(:,:,:)) 
     256         DO_3D( 0, 0, 0, 0, 1, jpk ) 
     257            zw3d(ji,jj,jk) = vn_25h(ji,jj,jk)*vmask(ji,jj,jk) + zmdi*(1.0-vmask(ji,jj,jk)) 
     258         END_3D 
    228259         CALL iom_put("vomecrty25h", zw3d  )   ! j-current 
    229          zw3d(:,:,:) = wn_25h(:,:,:)*wmask(:,:,:) + zmdi*(1.0-tmask(:,:,:)) 
     260         DO_3D( 0, 0, 0, 0, 1, jpk ) 
     261            zw3d(ji,jj,jk) = wn_25h(ji,jj,jk)*wmask(ji,jj,jk) + zmdi*(1.0-tmask(ji,jj,jk)) 
     262         END_3D 
    230263         CALL iom_put("vovecrtz25h", zw3d )   ! k-current 
    231264         ! Write vertical physics 
    232          zw3d(:,:,:) = avt_25h(:,:,:)*wmask(:,:,:) + zmdi*(1.0-tmask(:,:,:)) 
     265         DO_3D( 0, 0, 0, 0, 1, jpk ) 
     266            zw3d(ji,jj,jk) = avt_25h(ji,jj,jk)*wmask(ji,jj,jk) + zmdi*(1.0-tmask(ji,jj,jk)) 
     267         END_3D 
    233268         CALL iom_put("avt25h", zw3d )   ! diffusivity 
    234          zw3d(:,:,:) = avm_25h(:,:,:)*wmask(:,:,:) + zmdi*(1.0-tmask(:,:,:)) 
     269         DO_3D( 0, 0, 0, 0, 1, jpk ) 
     270            zw3d(ji,jj,jk) = avm_25h(ji,jj,jk)*wmask(ji,jj,jk) + zmdi*(1.0-tmask(ji,jj,jk)) 
     271         END_3D 
    235272         CALL iom_put("avm25h", zw3d)   ! viscosity 
    236273         IF( ln_zdftke ) THEN 
    237             zw3d(:,:,:) = en_25h(:,:,:)*wmask(:,:,:) + zmdi*(1.0-tmask(:,:,:)) 
     274            DO_3D( 0, 0, 0, 0, 1, jpk ) 
     275               zw3d(ji,jj,jk) = en_25h(ji,jj,jk)*wmask(ji,jj,jk) + zmdi*(1.0-tmask(ji,jj,jk)) 
     276            END_3D 
    238277            CALL iom_put("tke25h", zw3d)   ! tke 
    239278         ENDIF 
    240279         IF( ln_zdfgls ) THEN 
    241             zw3d(:,:,:) = en_25h(:,:,:)*wmask(:,:,:) + zmdi*(1.0-tmask(:,:,:)) 
     280            DO_3D( 0, 0, 0, 0, 1, jpk ) 
     281               zw3d(ji,jj,jk) = en_25h(ji,jj,jk)*wmask(ji,jj,jk) + zmdi*(1.0-tmask(ji,jj,jk)) 
     282            END_3D 
    242283            CALL iom_put("tke25h", zw3d)   ! tke 
    243             zw3d(:,:,:) = rmxln_25h(:,:,:)*wmask(:,:,:) + zmdi*(1.0-tmask(:,:,:)) 
     284            DO_3D( 0, 0, 0, 0, 1, jpk ) 
     285               zw3d(ji,jj,jk) = rmxln_25h(ji,jj,jk)*wmask(ji,jj,jk) + zmdi*(1.0-tmask(ji,jj,jk)) 
     286            END_3D 
    244287            CALL iom_put( "mxln25h",zw3d) 
    245288         ENDIF 
    246289         ! 
    247290         ! After the write reset the values to cnt=1 and sum values equal current value  
    248          tn_25h  (:,:,:) = ts (:,:,:,jp_tem,Kmm) 
    249          sn_25h  (:,:,:) = ts (:,:,:,jp_sal,Kmm) 
    250          sshn_25h(:,:)   = ssh(:,:,Kmm) 
    251          un_25h  (:,:,:) = uu  (:,:,:,Kmm) 
    252          vn_25h  (:,:,:) = vv  (:,:,:,Kmm) 
    253          wn_25h  (:,:,:) = ww  (:,:,:) 
    254          avt_25h (:,:,:) = avt (:,:,:) 
    255          avm_25h (:,:,:) = avm (:,:,:) 
     291         DO_3D( 0, 0, 0, 0, 1, jpk ) 
     292            tn_25h  (ji,jj,jk) = ts (ji,jj,jk,jp_tem,Kmm) 
     293            sn_25h  (ji,jj,jk) = ts (ji,jj,jk,jp_sal,Kmm) 
     294            un_25h  (ji,jj,jk) = uu (ji,jj,jk,Kmm) 
     295            vn_25h  (ji,jj,jk) = vv (ji,jj,jk,Kmm) 
     296            wn_25h  (ji,jj,jk) = ww (ji,jj,jk) 
     297            avt_25h (ji,jj,jk) = avt(ji,jj,jk) 
     298            avm_25h (ji,jj,jk) = avm(ji,jj,jk) 
     299         END_3D 
     300         DO_2D( 0, 0, 0, 0 ) 
     301            sshn_25h(ji,jj)    = ssh(ji,jj,Kmm) 
     302         END_2D 
    256303         IF( ln_zdftke ) THEN 
    257             en_25h(:,:,:) = en(:,:,:) 
     304            DO_3D( 0, 0, 0, 0, 1, jpk ) 
     305               en_25h(ji,jj,jk) = en(ji,jj,jk) 
     306            END_3D 
    258307         ENDIF 
    259308         IF( ln_zdfgls ) THEN 
    260             en_25h   (:,:,:) = en    (:,:,:) 
    261             rmxln_25h(:,:,:) = hmxl_n(:,:,:) 
     309            DO_3D( 0, 0, 0, 0, 1, jpk ) 
     310               en_25h   (ji,jj,jk) = en    (ji,jj,jk) 
     311               rmxln_25h(ji,jj,jk) = hmxl_n(ji,jj,jk) 
     312            END_3D 
    262313         ENDIF 
    263314         cnt_25h = 1 
  • NEMO/trunk/src/OCE/ZDF/zdfmxl.F90

    r15182 r15249  
    2626   PRIVATE 
    2727 
    28    PUBLIC   zdf_mxl, zdf_mxl_turb   ! called by zdfphy.F90 
     28   PUBLIC   zdf_mxl, zdf_mxl_turb, zdf_mxl_alloc   ! called by zdfphy.F90 
    2929 
    3030   INTEGER , PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) ::   nmln    !: number of level in the mixed layer (used by LDF, ZDF, TRD, TOP) 
     
    8686            IF(lwp) WRITE(numout,*) 'zdf_mxl : mixed layer depth' 
    8787            IF(lwp) WRITE(numout,*) '~~~~~~~ ' 
    88             !                             ! allocate zdfmxl arrays 
    89             IF( zdf_mxl_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_mxl : unable to allocate arrays' ) 
    9088         ENDIF 
    9189      ENDIF 
     
    111109      ! 
    112110      IF( .NOT.l_offline .AND. iom_use("mldr10_1") ) THEN 
    113          IF( ln_isfcav ) THEN  ;  CALL iom_put( "mldr10_1", hmlp - risfdep)   ! mixed layer thickness 
    114          ELSE                  ;  CALL iom_put( "mldr10_1", hmlp )            ! mixed layer depth 
    115          END IF 
     111         IF( .NOT. l_istiled .OR. ntile == nijtile ) THEN         ! Do only on the last tile 
     112            IF( ln_isfcav ) THEN  ;  CALL iom_put( "mldr10_1", hmlp - risfdep)   ! mixed layer thickness 
     113            ELSE                  ;  CALL iom_put( "mldr10_1", hmlp )            ! mixed layer depth 
     114            END IF 
     115         ENDIF 
    116116      ENDIF 
    117117      ! 
     
    154154      ! 
    155155      IF( .NOT.l_offline .AND. iom_use("mldkz5") ) THEN 
    156          IF( ln_isfcav ) THEN  ;  CALL iom_put( "mldkz5"  , hmld - risfdep )   ! turbocline thickness 
    157          ELSE                  ;  CALL iom_put( "mldkz5"  , hmld )             ! turbocline depth 
    158          END IF 
     156         IF( .NOT. l_istiled .OR. ntile == nijtile ) THEN         ! Do only on the last tile 
     157            IF( ln_isfcav ) THEN  ;  CALL iom_put( "mldkz5"  , hmld - risfdep )   ! turbocline thickness 
     158            ELSE                  ;  CALL iom_put( "mldkz5"  , hmld )             ! turbocline depth 
     159            END IF 
     160         ENDIF 
    159161      ENDIF 
    160162      ! 
  • NEMO/trunk/src/OCE/ZDF/zdfphy.F90

    r15182 r15249  
    146146         wi    (:,:,:) = 0._wp 
    147147      ENDIF 
     148      !                                  ! Initialise zdf_mxl arrays (only hmld as not set everywhere when nn_hls > 1) 
     149      IF( zdf_mxl_alloc() /= 0 )   CALL ctl_stop( 'STOP', 'zdf_mxl : unable to allocate arrays' ) 
     150      hmld(:,:) = 0._wp 
    148151      !                          !==  Background eddy viscosity and diffusivity  ==! 
    149152      IF( nn_avb == 0 ) THEN             ! Define avmb, avtb from namelist parameter 
  • NEMO/trunk/src/SWE/nemogcm.F90

    r15136 r15249  
    2525   USE isf_oce , ONLY : ln_isf           ! ice shelf 
    2626   USE trd_oce , ONLY : l_trddyn         ! dynamical trend logical 
     27   USE dia25h  , ONLY : ln_dia25h        ! 25h mean output 
    2728#if defined key_RK3 
    2829   USE stprk3         ! NEMO time-stepping               (stp_RK3   routine) 
     
    338339      !                                      ! Diagnostics 
    339340      IF( ln_diacfl    )   CALL dia_cfl_init    ! Initialise CFL diagnostics 
    340       !                                         ! Trends diag: switched off 
    341                            l_trddyn = .FALSE.        ! No trend diagnostics 
     341 
     342                           l_trddyn  = .FALSE.  ! No trend diagnostics 
     343                           ln_dia25h = .FALSE.  ! No 25h mean diagnostics (zdf_phy not used)- used in diawri 
    342344 
    343345      IF(lwp) WRITE(numout,cform_aaa)           ! Flag AAAAAAA 
  • NEMO/trunk/src/SWE/stpmlf.F90

    r14433 r15249  
    199199      CALL lbc_lnk( 'stp_MLF', uu(:,:,:,Nnn), 'U', -1., vv(:,:,:,Nnn), 'V', -1.,   &   !* local domain boundaries 
    200200         &                     uu(:,:,:,Naa), 'U', -1., vv(:,:,:,Naa), 'V', -1.    )      
     201      IF (nn_hls==2) CALL lbc_lnk( 'stp_MLF', r3u(:,:,Naa), 'U', 1., r3v(:,:,Naa), 'V', 1.) 
    201202 
    202203      !>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>>> 
  • NEMO/trunk/src/SWE/stprk3.F90

    r14834 r15249  
    172172      ! 
    173173      CALL lbc_lnk( 'stp_RK3', uu(:,:,:,Naa), 'U', -1., vv(:,:,:,Naa), 'V', -1. ) 
    174       IF (nn_hls==2) CALL lbc_lnk( 'stp_MLF', r3u(:,:,Naa), 'U', 1., r3v(:,:,Naa), 'U', 1.) 
     174      IF (nn_hls==2) CALL lbc_lnk( 'stp_RK3', r3u(:,:,Naa), 'U', 1., r3v(:,:,Naa), 'V', 1.) 
    175175      ! 
    176176      !                                 !==  Swap time levels  ==! 
     
    238238      ! 
    239239      CALL lbc_lnk( 'stp_RK3', uu(:,:,:,Naa), 'U', -1., vv(:,:,:,Naa), 'V', -1. ) 
    240       IF (nn_hls==2) CALL lbc_lnk( 'stp_MLF', r3u(:,:,Naa), 'U', 1., r3v(:,:,Naa), 'U', 1.) 
     240      IF (nn_hls==2) CALL lbc_lnk( 'stp_RK3', r3u(:,:,Naa), 'U', 1., r3v(:,:,Naa), 'V', 1.) 
    241241      ! 
    242242      !                                 !==  Swap time levels  ==! 
     
    302302      ! 
    303303      CALL lbc_lnk( 'stp_RK3', uu(:,:,:,Naa), 'U', -1., vv(:,:,:,Naa), 'V', -1. ) 
    304       IF (nn_hls==2) CALL lbc_lnk( 'stp_MLF', r3u(:,:,Naa), 'U', 1., r3v(:,:,Naa), 'U', 1.) 
     304      IF (nn_hls==2) CALL lbc_lnk( 'stp_RK3', r3u(:,:,Naa), 'U', 1., r3v(:,:,Naa), 'V', 1.) 
    305305      ! 
    306306      !                                 !==  Swap time levels  ==! 
Note: See TracChangeset for help on using the changeset viewer.