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

Changeset 15554 for NEMO/branches/2021


Ignore:
Timestamp:
2021-11-29T13:42:12+01:00 (2 years ago)
Author:
gsamson
Message:

adapt ABL routines to new potential temperature computation and reintroduce rn_vfac for ABL; ticket #2632

Location:
NEMO/branches/2021/ticket2632_r14588_theta_sbcblk
Files:
4 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2021/ticket2632_r14588_theta_sbcblk/cfgs/SHARED/namelist_ref

    r15548 r15554  
    331331   rn_ltra_min   =  4.5       ! tracers  nudging magnitude inside the ABL [hour] (~3 rn_Dt) 
    332332   rn_ltra_max   =  1.5       ! tracers  nudging magnitude above  the ABL [hour] (~1 rn_Dt) 
     333   rn_vfac       =  1.        ! multiplicative factor for ocean/ice velocity 
    333334   nn_amxl       =  0         ! mixing length: = 0 Deardorff 80 length-scale 
    334335                              !                = 1 length-scale based on the distance to the PBL height 
  • NEMO/branches/2021/ticket2632_r14588_theta_sbcblk/src/ABL/ablmod.F90

    r15548 r15554  
    8383      REAL(wp) , INTENT(inout), DIMENSION(:,:  ) ::   psen       ! Ch x Du 
    8484      REAL(wp) , INTENT(inout), DIMENSION(:,:  ) ::   pevp       ! Ce x Du 
    85       REAL(wp) , INTENT(inout), DIMENSION(:,:  ) ::   pwndm      ! ||uwnd|| 
     85      REAL(wp) , INTENT(inout), DIMENSION(:,:  ) ::   pwndm      ! ||uwnd - uoce|| 
    8686      REAL(wp) , INTENT(  out), DIMENSION(:,:  ) ::   plat       ! latent heat flux 
    8787      REAL(wp) , INTENT(  out), DIMENSION(:,:  ) ::   ptaui      ! taux 
     
    346346         DO jk = 3, jpkam1 
    347347            DO ji = 1, jpi 
    348                z_elem_a( ji, jk ) = - rDt_abl * Avm_abl( ji, jj, jk-1 ) / e3w_abl( jk-1 )  ! lower-diagonal 
    349                z_elem_c( ji, jk ) = - rDt_abl * Avm_abl( ji, jj, jk   ) / e3w_abl( jk   )  ! upper-diagonal 
    350                z_elem_b( ji, jk ) = e3t_abl(jk) - z_elem_a( ji, jk ) - z_elem_c( ji, jk )  !       diagonal 
    351             END DO 
    352          END DO 
    353  
    354          DO ji = 2, jpi   ! boundary conditions   (Avm_abl and pcd_du must be available at ji=jpi) 
     348               z_elem_a( ji, jk ) = - rDt_abl * Avm_abl( ji, jj, jk-1 ) / e3w_abl( jk-1 )   ! lower-diagonal 
     349               z_elem_c( ji, jk ) = - rDt_abl * Avm_abl( ji, jj, jk   ) / e3w_abl( jk   )   ! upper-diagonal 
     350               z_elem_b( ji, jk ) = e3t_abl(jk) - z_elem_a( ji, jk ) - z_elem_c( ji, jk )   !       diagonal 
     351            END DO 
     352         END DO 
     353 
     354         DO ji = 2, jpi   ! boundary conditions (Avm_abl and pcd_du must be available at ji=jpi) 
    355355            !++ Surface boundary condition 
    356356            z_elem_a( ji, 2 ) = 0._wp 
    357357            z_elem_c( ji, 2 ) = - rDt_abl * Avm_abl( ji, jj, 2 ) / e3w_abl( 2 ) 
    358358            ! 
    359             zztmp1  = pcd_du(ji, jj) 
    360             zztmp2  = 0.5_wp * pcd_du(ji, jj) * ( pssu(ji-1, jj) + pssu(ji,jj) ) 
     359            zztmp1 = pcd_du(ji, jj) 
     360            zztmp2 = 0.5_wp * pcd_du(ji, jj) * ( pssu(ji-1, jj) + pssu(ji, jj  ) ) * rn_vfac 
    361361#if defined key_si3 
    362362            zztmp1 = zztmp1 * pfrac_oce(ji,jj) + (1._wp - pfrac_oce(ji,jj)) * pcd_du_ice(ji, jj) 
    363             zzice  = 0.5_wp * ( pssu_ice(ji-1, jj) + pssu_ice(ji, jj) ) 
     363            zzice  = 0.5_wp * ( pssu_ice(ji-1, jj) + pssu_ice(ji, jj  ) ) * rn_vfac 
    364364            zztmp2 = zztmp2 * pfrac_oce(ji,jj) + (1._wp - pfrac_oce(ji,jj)) * pcd_du_ice(ji, jj) * zzice 
    365365#endif 
    366366            z_elem_b( ji,     2       ) = e3t_abl( 2 ) - z_elem_c( ji, 2 ) + rDt_abl * zztmp1 
    367             u_abl( ji, jj,    2, nt_a ) =      u_abl( ji, jj,    2, nt_a ) + rDt_abl * zztmp2 
     367            u_abl( ji, jj,    2, nt_a ) =         u_abl( ji, jj, 2, nt_a ) + rDt_abl * zztmp2 
    368368 
    369369            ! idealized test cases only 
     
    386386         !! Matrix inversion 
    387387         !! ---------------------------------------------------------- 
    388          !DO ji = 2, jpi 
    389          DO ji = 1, jpi  !!GS: TBI 
     388         DO ji = 1, jpi  !DO ji = 2, jpi !!GS: TBI 
    390389            zcff                 =   1._wp / z_elem_b( ji, 2 ) 
    391             zCF   (ji,   2     ) =  - zcff * z_elem_c( ji, 2 ) 
    392             u_abl (ji,jj,2,nt_a) =    zcff * u_abl(ji,jj,2,nt_a) 
     390            zCF   (ji,   2     ) =  - zcff * z_elem_c( ji,     2      ) 
     391            u_abl (ji,jj,2,nt_a) =    zcff * u_abl   ( ji, jj, 2, nt_a ) 
    393392         END DO 
    394393 
     
    420419         DO jk = 3, jpkam1 
    421420            DO ji = 1, jpi 
    422                z_elem_a( ji, jk ) = -rDt_abl * Avm_abl( ji, jj, jk-1 ) / e3w_abl( jk-1 )   ! lower-diagonal 
    423                z_elem_c( ji, jk ) = -rDt_abl * Avm_abl( ji, jj, jk   ) / e3w_abl( jk   )   ! upper-diagonal 
    424                z_elem_b( ji, jk ) = e3t_abl(jk) - z_elem_a( ji, jk ) - z_elem_c( ji, jk )                              !       diagonal 
     421               z_elem_a( ji, jk ) = - rDt_abl * Avm_abl( ji, jj, jk-1 ) / e3w_abl( jk-1 )   ! lower-diagonal 
     422               z_elem_c( ji, jk ) = - rDt_abl * Avm_abl( ji, jj, jk   ) / e3w_abl( jk   )   ! upper-diagonal 
     423               z_elem_b( ji, jk ) = e3t_abl(jk) - z_elem_a( ji, jk ) - z_elem_c( ji, jk )   !       diagonal 
    425424            END DO 
    426425         END DO 
     
    432431            ! 
    433432            zztmp1 = pcd_du(ji, jj) 
    434             zztmp2 = 0.5_wp * pcd_du(ji, jj) * ( pssv(ji, jj) + pssv(ji, jj-1) ) 
     433            zztmp2 = 0.5_wp * pcd_du(ji, jj) * ( pssv(ji, jj) + pssv(ji, jj-1) ) * rn_vfac 
    435434#if defined key_si3 
    436435            zztmp1 = zztmp1 * pfrac_oce(ji,jj) + (1._wp - pfrac_oce(ji,jj)) * pcd_du_ice(ji, jj) 
    437             zzice  = 0.5_wp * ( pssv_ice(ji, jj) + pssv_ice(ji, jj-1) ) 
     436            zzice  = 0.5_wp * ( pssv_ice(ji, jj) + pssv_ice(ji, jj-1) ) * rn_vfac 
    438437            zztmp2 = zztmp2 * pfrac_oce(ji,jj) + (1._wp - pfrac_oce(ji,jj)) * pcd_du_ice(ji, jj) * zzice 
    439438#endif 
     
    461460         !! ---------------------------------------------------------- 
    462461         DO ji = 1, jpi 
    463             zcff                 =  1._wp / z_elem_b( ji, 2 ) 
    464             zCF   (ji,   2     ) =   - zcff * z_elem_c( ji,     2       ) 
    465             v_abl (ji,jj,2,nt_a) =     zcff * v_abl   ( ji, jj, 2, nt_a ) 
     462            zcff                 =   1._wp / z_elem_b( ji, 2 ) 
     463            zCF   (ji,   2     ) =  - zcff * z_elem_c( ji,     2       ) 
     464            v_abl (ji,jj,2,nt_a) =    zcff * v_abl   ( ji, jj, 2, nt_a ) 
    466465         END DO 
    467466 
     
    501500                  &  + jp_alp1_dyn * zsig    + jp_alp0_dyn 
    502501               zcff  = (1._wp-zmsk) + zmsk * zcff2 * rDt_abl   ! zcff = 1 for masked points 
     502                                                               ! rn_Dt = rDt_abl / nn_fsbc 
    503503               zcff  = zcff * rest_eq(ji,jj) 
    504504               u_abl( ji, jj, jk, nt_a ) = (1._wp - zcff ) *  u_abl( ji, jj, jk, nt_a )   & 
     
    523523               &  + jp_alp1_tra * zsig    + jp_alp0_tra 
    524524            zcff  = (1._wp-zmsk) + zmsk * zcff2 * rDt_abl   ! zcff = 1 for masked points 
     525                                                            ! rn_Dt = rDt_abl / nn_fsbc 
    525526            tq_abl( ji, jj, jk, nt_a, jp_ta ) = (1._wp - zcff ) * tq_abl( ji, jj, jk, nt_a, jp_ta )   & 
    526527               &                                       + zcff   * pt_dta( ji, jj, jk ) 
     
    538539      ! 
    539540      CALL lbc_lnk( 'ablmod',  u_abl(:,:,:,nt_a      ), 'T', -1._wp,  v_abl(:,:,:,nt_a)      , 'T', -1._wp                            ) 
    540       CALL lbc_lnk( 'ablmod', tq_abl(:,:,:,nt_a,jp_ta), 'T', 1._wp , tq_abl(:,:,:,nt_a,jp_qa), 'T',  1._wp , kfillmode = jpfillnothing )   ! ++++ this should not be needed... 
     541      !CALL lbc_lnk( 'ablmod', tq_abl(:,:,:,nt_a,jp_ta), 'T',  1._wp, tq_abl(:,:,:,nt_a,jp_qa), 'T',  1._wp, kfillmode = jpfillnothing )   ! ++++ this should not be needed... 
    541542      ! 
    542543#if defined key_xios 
     
    557558      END IF 
    558559      IF( ln_hpgls_frc ) THEN 
    559          IF ( iom_use("uz1_geo") ) CALL iom_put ( "uz1_geo",  pgu_dta(:,:,2)/MAX(fft_abl(:,:),2.5e-5_wp) ) 
    560          IF ( iom_use("vz1_geo") ) CALL iom_put ( "vz1_geo", -pgv_dta(:,:,2)/MAX(fft_abl(:,:),2.5e-5_wp) ) 
     560         IF ( iom_use("uz1_geo") ) CALL iom_put ( "uz1_geo", -pgv_dta(:,:,2)/MAX( ABS( fft_abl(:,:)), 2.5e-5_wp)*fft_abl(:,:)/ABS(fft_abl(:,:)) ) 
     561         IF ( iom_use("vz1_geo") ) CALL iom_put ( "vz1_geo",  pgu_dta(:,:,2)/MAX( ABS( fft_abl(:,:)), 2.5e-5_wp)*fft_abl(:,:)/ABS(fft_abl(:,:)) ) 
    561562      END IF 
    562563      ! 3D (all ABL levels) 
     
    580581      END IF 
    581582      IF( ln_hpgls_frc ) THEN 
    582          IF ( iom_use("u_geo") ) CALL iom_put ( "u_geo",  pgu_dta(:,:,2:jpka)/MAX( RESHAPE( fft_abl(:,:), (/jpi,jpj,jpka-1/), fft_abl(:,:)), 2.5e-5_wp) ) 
    583          IF ( iom_use("v_geo") ) CALL iom_put ( "v_geo", -pgv_dta(:,:,2:jpka)/MAX( RESHAPE( fft_abl(:,:), (/jpi,jpj,jpka-1/), fft_abl(:,:)), 2.5e-5_wp) ) 
     583         IF ( iom_use("u_geo") ) CALL iom_put ( "u_geo", -pgv_dta(:,:,2:jpka)/MAX( ABS( RESHAPE( fft_abl(:,:), (/jpi,jpj,jpka-1/), fft_abl(:,:))), 2.5e-5_wp) * RESHAPE( fft_abl(:,:)/ABS(fft_abl(:,:)), (/jpi,jpj,jpka-1/), fft_abl(:,:)/ABS(fft_abl(:,:))) ) 
     584         IF ( iom_use("v_geo") ) CALL iom_put ( "v_geo",  pgu_dta(:,:,2:jpka)/MAX( ABS( RESHAPE( fft_abl(:,:), (/jpi,jpj,jpka-1/), fft_abl(:,:))), 2.5e-5_wp) * RESHAPE( fft_abl(:,:)/ABS(fft_abl(:,:)), (/jpi,jpj,jpka-1/), fft_abl(:,:)/ABS(fft_abl(:,:))) ) 
    584585      END IF 
    585586#endif 
     
    598599         rhoa( ji, jj ) = zcff 
    599600      END_2D 
    600       !!GS: debug only; to be removed 
    601       CALL iom_put ( "pres_zu", zpre(:,:) )  
    602       CALL iom_put ( "tabs_zu", ztabs(:,:) )  
    603601 
    604602      DO_2D( nn_hls-1, nn_hls, nn_hls-1, nn_hls ) 
    605          zwnd_i(ji,jj) = u_abl(ji  ,jj,2,nt_a) - 0.5_wp * ( pssu(ji  ,jj) + pssu(ji-1,jj) ) 
    606          zwnd_j(ji,jj) = v_abl(ji,jj  ,2,nt_a) - 0.5_wp * ( pssv(ji,jj  ) + pssv(ji,jj-1) ) 
     603         zwnd_i(ji,jj) = u_abl(ji  ,jj,2,nt_a) - 0.5_wp * ( pssu(ji  ,jj) + pssu(ji-1,jj) ) * rn_vfac 
     604         zwnd_j(ji,jj) = v_abl(ji,jj  ,2,nt_a) - 0.5_wp * ( pssv(ji,jj  ) + pssv(ji,jj-1) ) * rn_vfac 
    607605      END_2D 
    608606      ! 
     
    667665         ptaui_ice(ji,jj) = 0.5_wp * (  rhoa(ji+1,jj) * pCd_du_ice(ji+1,jj)             & 
    668666            &                      +    rhoa(ji  ,jj) * pCd_du_ice(ji  ,jj)  )          & 
    669             &         * ( zztmp1 - pssu_ice(ji,jj) ) 
     667            &         * ( zztmp1 - pssu_ice(ji,jj) * rn_vfac ) 
    670668         ptauj_ice(ji,jj) = 0.5_wp * (  rhoa(ji,jj+1) * pCd_du_ice(ji,jj+1)             & 
    671669            &                      +    rhoa(ji,jj  ) * pCd_du_ice(ji,jj  )  )          & 
    672             &         * ( zztmp2 - pssv_ice(ji,jj) ) 
     670            &         * ( zztmp2 - pssv_ice(ji,jj) * rn_vfac ) 
    673671      END_2D 
    674672      CALL lbc_lnk( 'ablmod', ptaui_ice, 'U', -1.0_wp, ptauj_ice,'V', -1.0_wp ) 
  • NEMO/branches/2021/ticket2632_r14588_theta_sbcblk/src/ABL/par_abl.F90

    r13214 r15554  
    6161   REAL(wp), PUBLIC            ::   rn_ltra_max                   !: namelist parameter 
    6262   REAL(wp), PUBLIC            ::   rn_Ric                        !: critical Richardson number  
     63   REAL(wp), PUBLIC            ::   rn_vfac                       !: multiplicative factor for ocean/ice velocity 
    6364 
    6465   !!--------------------------------------------------------------------- 
  • NEMO/branches/2021/ticket2632_r14588_theta_sbcblk/src/ABL/sbcabl.F90

    r14592 r15554  
    7272         &                 rn_ldyn_min , rn_ldyn_max, rn_ltra_min, rn_ltra_max,   & 
    7373         &                 nn_amxl, rn_Cm, rn_Ct, rn_Ce, rn_Ceps, rn_Rod, rn_Ric, & 
    74          &                 ln_smth_pblh 
     74         &                 rn_vfac, ln_smth_pblh 
    7575      !!--------------------------------------------------------------------- 
    7676 
     
    370370#endif 
    371371            &                                                                  ) 
     372 
    372373         !!------------------------------------------------------------------------------------------- 
    373374         !! 4 - Finalize flux computation using ABL variables at (n+1), nt_n corresponds to (n+1) since 
Note: See TracChangeset for help on using the changeset viewer.