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 8898 for branches/2017/dev_r8624_AGRIF3_VVL/NEMOGCM/NEMO/OPA_SRC/DYN – NEMO

Ignore:
Timestamp:
2017-12-05T14:06:56+01:00 (6 years ago)
Author:
jchanut
Message:

AGRIF+VVL: Merge with free surface changes - #1965

Location:
branches/2017/dev_r8624_AGRIF3_VVL/NEMOGCM/NEMO/OPA_SRC/DYN
Files:
4 edited

Legend:

Unmodified
Added
Removed
  • branches/2017/dev_r8624_AGRIF3_VVL/NEMOGCM/NEMO/OPA_SRC/DYN/dynnxt.F90

    r8741 r8898  
    209209            ! (used as a now filtered scale factor until the swap) 
    210210            ! ---------------------------------------------------- 
    211             IF( ln_dynspg_ts .AND. ln_bt_fw ) THEN    ! No asselin filtering on thicknesses if forward time splitting 
    212                e3t_b(:,:,1:jpkm1) = e3t_n(:,:,1:jpkm1) 
    213             ELSE 
    214                DO jk = 1, jpkm1 
    215                   e3t_b(:,:,jk) = e3t_n(:,:,jk) + atfp * ( e3t_b(:,:,jk) - 2._wp * e3t_n(:,:,jk) + e3t_a(:,:,jk) ) 
     211            DO jk = 1, jpkm1 
     212               e3t_b(:,:,jk) = e3t_n(:,:,jk) + atfp * ( e3t_b(:,:,jk) - 2._wp * e3t_n(:,:,jk) + e3t_a(:,:,jk) ) 
     213            END DO 
     214            ! Add volume filter correction: compatibility with tracer advection scheme 
     215            ! => time filter + conservation correction (only at the first level) 
     216            zcoef = atfp * rdt * r1_rau0 
     217            IF ( .NOT. ln_isf ) THEN   ! if no ice shelf melting 
     218               e3t_b(:,:,1) = e3t_b(:,:,1) - zcoef * ( emp_b(:,:) - emp(:,:) & 
     219                              &                      - rnf_b(:,:) + rnf(:,:) ) * tmask(:,:,1) 
     220            ELSE                     ! if ice shelf melting 
     221               DO jj = 1, jpj 
     222                  DO ji = 1, jpi 
     223                     ikt = mikt(ji,jj) 
     224                     e3t_b(ji,jj,ikt) = e3t_b(ji,jj,ikt) - zcoef * (  emp_b   (ji,jj) - emp   (ji,jj)  & 
     225                        &                                           - rnf_b   (ji,jj) + rnf   (ji,jj)  & 
     226                        &                                           + fwfisf_b(ji,jj) - fwfisf(ji,jj)  ) * tmask(ji,jj,ikt) 
     227                  END DO 
    216228               END DO 
    217                ! Add volume filter correction: compatibility with tracer advection scheme 
    218                ! => time filter + conservation correction (only at the first level) 
    219                zcoef = atfp * rdt * r1_rau0 
    220                IF ( .NOT. ln_isf ) THEN   ! if no ice shelf melting 
    221                   e3t_b(:,:,1) = e3t_b(:,:,1) - zcoef * ( emp_b(:,:) - emp(:,:) & 
    222                                  &                      - rnf_b(:,:) + rnf(:,:) ) * tmask(:,:,1) 
    223                ELSE                     ! if ice shelf melting 
    224                   DO jj = 1, jpj 
    225                      DO ji = 1, jpi 
    226                         ikt = mikt(ji,jj) 
    227                         e3t_b(ji,jj,ikt) = e3t_b(ji,jj,ikt) - zcoef * (  emp_b   (ji,jj) - emp   (ji,jj)  & 
    228                            &                                           - rnf_b   (ji,jj) + rnf   (ji,jj)  & 
    229                            &                                           + fwfisf_b(ji,jj) - fwfisf(ji,jj)  ) * tmask(ji,jj,ikt) 
    230                      END DO 
    231                   END DO 
    232                END IF 
    233             ENDIF 
     229            END IF 
    234230            ! 
    235231            IF( ln_dynadv_vec ) THEN      ! Asselin filter applied on velocity 
  • branches/2017/dev_r8624_AGRIF3_VVL/NEMOGCM/NEMO/OPA_SRC/DYN/dynspg.F90

    r7753 r8898  
    183183      NAMELIST/namdyn_spg/ ln_dynspg_exp       , ln_dynspg_ts,   & 
    184184      &                    ln_bt_fw, ln_bt_av  , ln_bt_auto  ,   & 
    185       &                    nn_baro , rn_bt_cmax, nn_bt_flt 
     185      &                    nn_baro , rn_bt_cmax, nn_bt_flt, rn_bt_alpha 
    186186      !!---------------------------------------------------------------------- 
    187187      ! 
  • branches/2017/dev_r8624_AGRIF3_VVL/NEMOGCM/NEMO/OPA_SRC/DYN/dynspg_ts.F90

    r8762 r8898  
    6767   PUBLIC ts_rst            !    "      "     "    " 
    6868 
    69    INTEGER, SAVE :: icycle  ! Number of barotropic sub-steps for each internal step nn_baro <= 2.5 nn_baro 
    70    REAL(wp),SAVE :: rdtbt   ! Barotropic time step 
     69   INTEGER, SAVE :: icycle      ! Number of barotropic sub-steps for each internal step nn_baro <= 2.5 nn_baro 
     70   REAL(wp),SAVE :: rdtbt       ! Barotropic time step 
    7171 
    7272   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:) ::   wgtbtp1, wgtbtp2   !: 1st & 2nd weights used in time filtering of barotropic fields 
     
    151151      REAL(wp) ::   zhura, zhvra          !   -      - 
    152152      REAL(wp) ::   za0, za1, za2, za3    !   -      - 
    153       ! 
     153      REAL(wp) ::   zepsilon, zgamma            !   -      - 
    154154      REAL(wp), POINTER, DIMENSION(:,:) :: zsshp2_e 
    155155      REAL(wp), POINTER, DIMENSION(:,:) :: zu_trd, zv_trd, zu_frc, zv_frc, zssh_frc 
     
    759759           za3= 0._wp               
    760760         ELSE                               ! AB3-AM4 Coefficients; bet=0.281105 ; eps=0.013 ; gam=0.0880  
    761            za0=0.614_wp                     ! za0 = 1/2 +   gam + 2*eps     
    762            za1=0.285_wp                     ! za1 = 1/2 - 2*gam - 3*eps  
    763            za2=0.088_wp                     ! za2 = gam 
    764            za3=0.013_wp                     ! za3 = eps 
     761            IF (rn_bt_alpha==0._wp) THEN 
     762               za0=0.614_wp                 ! za0 = 1/2 +   gam + 2*eps 
     763               za1=0.285_wp                 ! za1 = 1/2 - 2*gam - 3*eps 
     764               za2=0.088_wp                 ! za2 = gam 
     765               za3=0.013_wp                 ! za3 = eps 
     766            ELSE 
     767               zepsilon = 0.00976186_wp - 0.13451357_wp * rn_bt_alpha 
     768               zgamma = 0.08344500_wp - 0.51358400_wp * rn_bt_alpha 
     769               za0 = 0.5_wp + zgamma + 2._wp * rn_bt_alpha + 2._wp * zepsilon 
     770               za1 = 1._wp - za0 - zgamma - zepsilon 
     771               za2 = zgamma 
     772               za3 = zepsilon 
     773            ENDIF  
    765774         ENDIF 
    766775         ! 
     
    10241033         zwy(:,:) = vn_adv(:,:) 
    10251034         IF( .NOT.( kt == nit000 .AND. neuler==0 ) ) THEN 
    1026             un_adv(:,:) = z1_2 * ( ub2_b(:,:) + zwx(:,:) )  
    1027             vn_adv(:,:) = z1_2 * ( vb2_b(:,:) + zwy(:,:) ) 
    1028          END IF 
     1035            un_adv(:,:) = z1_2 * ( ub2_b(:,:) + zwx(:,:) - atfp * un_bf(:,:) ) 
     1036            vn_adv(:,:) = z1_2 * ( vb2_b(:,:) + zwy(:,:) - atfp * vn_bf(:,:) ) 
     1037            ! 
     1038            ! Update corrective fluxes for next time step: 
     1039            un_bf(:,:)  = atfp * un_bf(:,:) + (zwx(:,:) - ub2_b(:,:)) 
     1040            vn_bf(:,:)  = atfp * vn_bf(:,:) + (zwy(:,:) - vb2_b(:,:)) 
     1041         ELSE 
     1042            un_bf(:,:) = 0._wp 
     1043            vn_bf(:,:) = 0._wp  
     1044         END IF          
    10291045         ! Save integrated transport for next computation 
    10301046         ub2_b(:,:) = zwx(:,:) 
     
    11981214         CALL iom_get( numror, jpdom_autoglo, 'ub2_b'  , ub2_b  (:,:) )    
    11991215         CALL iom_get( numror, jpdom_autoglo, 'vb2_b'  , vb2_b  (:,:) )  
     1216         CALL iom_get( numror, jpdom_autoglo, 'un_bf'  , un_bf  (:,:) )    
     1217         CALL iom_get( numror, jpdom_autoglo, 'vn_bf'  , vn_bf  (:,:) )  
    12001218         IF( .NOT.ln_bt_av ) THEN 
    12011219            CALL iom_get( numror, jpdom_autoglo, 'sshbb_e'  , sshbb_e(:,:) )    
     
    12171235         CALL iom_rstput( kt, nitrst, numrow, 'ub2_b'   , ub2_b  (:,:) ) 
    12181236         CALL iom_rstput( kt, nitrst, numrow, 'vb2_b'   , vb2_b  (:,:) ) 
     1237         CALL iom_rstput( kt, nitrst, numrow, 'un_bf'   , un_bf  (:,:) ) 
     1238         CALL iom_rstput( kt, nitrst, numrow, 'vn_bf'   , vn_bf  (:,:) ) 
    12191239         ! 
    12201240         IF (.NOT.ln_bt_av) THEN 
     
    12951315#if defined key_agrif 
    12961316      ! Restrict the use of Agrif to the forward case only 
    1297       IF( .NOT.ln_bt_fw .AND. .NOT.Agrif_Root() )   CALL ctl_stop( 'AGRIF not implemented if ln_bt_fw=.FALSE.' ) 
     1317!!!      IF( .NOT.ln_bt_fw .AND. .NOT.Agrif_Root() )   CALL ctl_stop( 'AGRIF not implemented if ln_bt_fw=.FALSE.' ) 
    12981318#endif 
    12991319      ! 
     
    13111331      IF(lwp) WRITE(numout,*) '     Maximum Courant number is   :', zcmax 
    13121332      ! 
     1333      IF(lwp) WRITE(numout,*)    '     Time diffusion parameter rn_bt_alpha: ', rn_bt_alpha 
     1334      IF ((ln_bt_av.AND.nn_bt_flt/=0).AND.(rn_bt_alpha>0._wp)) THEN 
     1335         CALL ctl_stop( 'dynspg_ts ERROR: if rn_bt_alpha > 0, remove temporal averaging' ) 
     1336      ENDIF 
     1337      ! 
    13131338      IF( .NOT.ln_bt_av .AND. .NOT.ln_bt_fw ) THEN 
    13141339         CALL ctl_stop( 'dynspg_ts ERROR: No time averaging => only forward integration is possible' ) 
  • branches/2017/dev_r8624_AGRIF3_VVL/NEMOGCM/NEMO/OPA_SRC/DYN/sshwzv.F90

    r8741 r8898  
    260260      ENDIF 
    261261      !              !==  Euler time-stepping: no filter, just swap  ==! 
    262       IF(  ( neuler == 0 .AND. kt == nit000 ) .OR.    & 
    263          & ( ln_bt_fw    .AND. ln_dynspg_ts )      ) THEN  
     262      IF  ( neuler == 0 .AND. kt == nit000 ) THEN 
    264263         sshb(:,:) = sshn(:,:)                              ! before <-- now 
    265264         sshn(:,:) = ssha(:,:)                              ! now    <-- after  (before already = now) 
Note: See TracChangeset for help on using the changeset viewer.