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

Changeset 13398


Ignore:
Timestamp:
2020-08-14T12:07:46+02:00 (4 years ago)
Author:
agn
Message:

include option of limiting hmle

Location:
NEMO/branches/2019/dev_r11078_OSMOSIS_IMMERSE_Nurser_4.0
Files:
2 edited

Legend:

Unmodified
Added
Removed
  • NEMO/branches/2019/dev_r11078_OSMOSIS_IMMERSE_Nurser_4.0/cfgs/SHARED/namelist_ref

    r12314 r13398  
    55!! namelists    2 - Surface boundary (namsbc, namsbc_flx, namsbc_blk, namsbc_cpl, 
    66!!                                    namsbc_sas, namtra_qsr, namsbc_rnf, 
    7 !!                                    namsbc_isf, namsbc_iscpl, namsbc_apr,  
     7!!                                    namsbc_isf, namsbc_iscpl, namsbc_apr, 
    88!!                                    namsbc_ssr, namsbc_wave, namberg) 
    99!!              3 - lateral boundary (namlbc, namagrif, nambdy, nambdy_tide) 
     
    8888      cn_domcfg = "domain_cfg"  ! domain configuration filename 
    8989      ! 
    90       ln_closea    = .false.    !  T => keep closed seas (defined by closea_mask field) in the   
     90      ln_closea    = .false.    !  T => keep closed seas (defined by closea_mask field) in the 
    9191      !                         !       domain and apply special treatment of freshwater fluxes. 
    92       !                         !  F => suppress closed seas (defined by closea_mask field)  
     92      !                         !  F => suppress closed seas (defined by closea_mask field) 
    9393      !                         !       from the bathymetry at runtime. 
    9494      !                         !  If closea_mask field doesn't exist in the domain_cfg file 
     
    106106   ln_tsd_init = .false.         !  ocean initialisation 
    107107   ln_tsd_dmp  = .false.         !  T-S restoring   (see namtra_dmp) 
    108     
     108 
    109109   cn_dir      = './'      !  root directory for the T-S data location 
    110110   !___________!_________________________!___________________!___________!_____________!________!___________!__________________!__________!_______________! 
     
    195195   nn_fsbc     = 2         !  frequency of SBC module call 
    196196      !                    !  (control sea-ice & iceberg model call) 
    197                      ! Type of air-sea fluxes  
     197                     ! Type of air-sea fluxes 
    198198   ln_usr      = .false.   !  user defined formulation                  (T => check usrdef_sbc) 
    199199   ln_flx      = .false.   !  flux formulation                          (T => fill namsbc_flx ) 
     
    205205      !                    !  =0 no opa-sas OASIS coupling: default single executable config. 
    206206      !                    !  =1 opa-sas OASIS coupling: multi executable config., OPA component 
    207       !                    !  =2 opa-sas OASIS coupling: multi executable config., SAS component  
     207      !                    !  =2 opa-sas OASIS coupling: multi executable config., SAS component 
    208208                     ! Sea-ice : 
    209    nn_ice      = 0         !  =0 no ice boundary condition     
     209   nn_ice      = 0         !  =0 no ice boundary condition 
    210210      !                    !  =1 use observed ice-cover                 (  => fill namsbc_iif ) 
    211211      !                    !  =2 or 3 automatically for SI3 or CICE    ("key_si3" or "key_cice") 
     
    213213   ln_ice_embd = .false.   !  =T embedded sea-ice (pressure + mass and salt exchanges) 
    214214      !                    !  =F levitating ice (no pressure, mass and salt exchanges) 
    215                      ! Misc. options of sbc :  
     215                     ! Misc. options of sbc : 
    216216   ln_traqsr   = .false.   !  Light penetration in the ocean            (T => fill namtra_qsr) 
    217217   ln_dm2dc    = .false.   !  daily mean to diurnal cycle on short wave 
     
    225225   ln_wave     = .false.   !  Activate coupling with wave  (T => fill namsbc_wave) 
    226226   ln_cdgw     = .false.   !  Neutral drag coefficient read from wave model (T => ln_wave=.true. & fill namsbc_wave) 
    227    ln_sdw      = .false.   !  Read 2D Surf Stokes Drift & Computation of 3D stokes drift (T => ln_wave=.true. & fill namsbc_wave)  
     227   ln_sdw      = .false.   !  Read 2D Surf Stokes Drift & Computation of 3D stokes drift (T => ln_wave=.true. & fill namsbc_wave) 
    228228   nn_sdrift   =  0        !  Parameterization for the calculation of 3D-Stokes drift from the surface Stokes drift 
    229229      !                    !   = 0 Breivik 2015 parameterization: v_z=v_0*[exp(2*k*z)/(1-8*k*z)] 
     
    379379   nn_chldta   =      0       !  RGB : Chl data (=1) or cst value (=0) 
    380380   rn_si1      =   23.0       !  2BD : longest depth of extinction 
    381     
     381 
    382382   cn_dir      = './'      !  root directory for the chlorophyl data location 
    383383   !___________!_________________________!___________________!___________!_____________!________!___________!__________________!__________!_______________! 
     
    443443/ 
    444444!----------------------------------------------------------------------- 
    445 &namsbc_isf    !  Top boundary layer (ISF)                              (ln_isfcav =T : read (ln_read_cfg=T)  
     445&namsbc_isf    !  Top boundary layer (ISF)                              (ln_isfcav =T : read (ln_read_cfg=T) 
    446446!-----------------------------------------------------------------------             or set or usr_def_zgr ) 
    447    !                 ! type of top boundary layer  
     447   !                 ! type of top boundary layer 
    448448   nn_isf      = 1         !  ice shelf melting/freezing 
    449                            !  1 = presence of ISF   ;  2 = bg03 parametrisation  
     449                           !  1 = presence of ISF   ;  2 = bg03 parametrisation 
    450450                           !  3 = rnf file for ISF  ;  4 = ISF specified freshwater flux 
    451451                           !  options 1 and 4 need ln_isfcav = .true. (domzgr) 
     
    470470!* nn_isf = 3 case 
    471471   sn_rnfisf   = 'rnfisf'    ,         -12.      ,'sofwfisf' ,  .false.    , .true.  , 'yearly'  ,    ''    ,   ''     ,    '' 
    472 !* nn_isf = 2 and 3 cases  
     472!* nn_isf = 2 and 3 cases 
    473473   sn_depmax_isf ='rnfisf'   ,         -12.      ,'sozisfmax',  .false.    , .true.  , 'yearly'  ,    ''    ,   ''     ,    '' 
    474474   sn_depmin_isf ='rnfisf'   ,         -12.      ,'sozisfmin',  .false.    , .true.  , 'yearly'  ,    ''    ,   ''     ,    '' 
     
    477477/ 
    478478!----------------------------------------------------------------------- 
    479 &namsbc_iscpl  !   land ice / ocean coupling option                     (ln_isfcav =T : read (ln_read_cfg=T)  
     479&namsbc_iscpl  !   land ice / ocean coupling option                     (ln_isfcav =T : read (ln_read_cfg=T) 
    480480!-----------------------------------------------------------------------             or set or usr_def_zgr ) 
    481481   nn_drown    = 10        ! number of iteration of the extrapolation loop (fill the new wet cells) 
     
    577577         ln_read_load  = .false.               ! Or read load potential from file 
    578578            cn_tide_load = 'tide_LOAD_grid_T.nc'  ! filename for load potential 
    579             !       
     579            ! 
    580580      ln_tide_ramp  = .false.               !  Use linear ramp for tides at startup 
    581581         rdttideramp   =    0.                 !  ramp duration in days 
     
    656656   filtide          = 'bdydta/amm12_bdytide_'   !  file name root of tidal forcing files 
    657657   ln_bdytide_2ddta = .false.                   ! 
    658    ln_bdytide_conj  = .false.                   !  
     658   ln_bdytide_conj  = .false.                   ! 
    659659/ 
    660660 
     
    683683!----------------------------------------------------------------------- 
    684684   rn_Cd0      =  1.e-3    !  drag coefficient [-] 
    685    rn_Uc0      =  0.4      !  ref. velocity [m/s] (linear drag=Cd0*Uc0)  
     685   rn_Uc0      =  0.4      !  ref. velocity [m/s] (linear drag=Cd0*Uc0) 
    686686   rn_Cdmax    =  0.1      !  drag value maximum [-] (logarithmic drag) 
    687687   rn_ke0      =  2.5e-3   !  background kinetic energy  [m2/s2] (non-linear cases) 
     
    694694!----------------------------------------------------------------------- 
    695695   rn_Cd0      =  1.e-3    !  drag coefficient [-] 
    696    rn_Uc0      =  0.4      !  ref. velocity [m/s] (linear drag=Cd0*Uc0)  
     696   rn_Uc0      =  0.4      !  ref. velocity [m/s] (linear drag=Cd0*Uc0) 
    697697   rn_Cdmax    =  0.1      !  drag value maximum [-] (logarithmic drag) 
    698698   rn_ke0      =  2.5e-3   !  background kinetic energy  [m2/s2] (non-linear cases) 
     
    761761      nn_cen_v   =  4            !  =2/4, vertical   2nd order CEN / 4th order COMPACT 
    762762   ln_traadv_fct = .false. !  FCT scheme 
    763       nn_fct_h   =  2            !  =2/4, horizontal 2nd / 4th order  
    764       nn_fct_v   =  2            !  =2/4, vertical   2nd / COMPACT 4th order  
     763      nn_fct_h   =  2            !  =2/4, horizontal 2nd / 4th order 
     764      nn_fct_v   =  2            !  =2/4, vertical   2nd / COMPACT 4th order 
    765765   ln_traadv_mus = .false. !  MUSCL scheme 
    766766      ln_mus_ups = .false.       !  use upstream scheme near river mouths 
     
    783783   ln_traldf_triad = .false.   !  iso-neutral (triad    operator) 
    784784   ! 
    785    !                       !  iso-neutral options:         
     785   !                       !  iso-neutral options: 
    786786   ln_traldf_msc   = .false.   !  Method of Stabilizing Correction      (both operators) 
    787787   rn_slpmax       =  0.01     !  slope limit                           (both operators) 
     
    793793   nn_aht_ijk_t    = 0         !  space/time variation of eddy coefficient: 
    794794      !                             !   =-20 (=-30)    read in eddy_diffusivity_2D.nc (..._3D.nc) file 
    795       !                             !   =  0           constant  
    796       !                             !   = 10 F(k)      =ldf_c1d  
    797       !                             !   = 20 F(i,j)    =ldf_c2d  
     795      !                             !   =  0           constant 
     796      !                             !   = 10 F(k)      =ldf_c1d 
     797      !                             !   = 20 F(i,j)    =ldf_c2d 
    798798      !                             !   = 21 F(i,j,t)  =Treguier et al. JPO 1997 formulation 
    799799      !                             !   = 30 F(i,j,k)  =ldf_c2d * ldf_c1d 
    800800      !                             !   = 31 F(i,j,k,t)=F(local velocity and grid-spacing) 
    801       !                        !  time invariant coefficients:  aht0 = 1/2  Ud*Ld   (lap case)  
     801      !                        !  time invariant coefficients:  aht0 = 1/2  Ud*Ld   (lap case) 
    802802      !                             !                           or   = 1/12 Ud*Ld^3 (blp case) 
    803803      rn_Ud        = 0.01           !  lateral diffusive velocity [m/s] (nn_aht_ijk_t= 0, 10, 20, 30) 
     
    825825      nn_aei_ijk_t    = 0           !  space/time variation of eddy coefficient: 
    826826      !                             !   =-20 (=-30)    read in eddy_induced_velocity_2D.nc (..._3D.nc) file 
    827       !                             !   =  0           constant  
    828       !                             !   = 10 F(k)      =ldf_c1d  
    829       !                             !   = 20 F(i,j)    =ldf_c2d  
     827      !                             !   =  0           constant 
     828      !                             !   = 10 F(k)      =ldf_c1d 
     829      !                             !   = 20 F(i,j)    =ldf_c2d 
    830830      !                             !   = 21 F(i,j,t)  =Treguier et al. JPO 1997 formulation 
    831831      !                             !   = 30 F(i,j,k)  =ldf_c2d * ldf_c1d 
    832       !                        !  time invariant coefficients:  aei0 = 1/2  Ue*Le  
     832      !                        !  time invariant coefficients:  aei0 = 1/2  Ue*Le 
    833833      rn_Ue        = 0.02           !  lateral diffusive velocity [m/s] (nn_aht_ijk_t= 0, 10, 20, 30) 
    834834      rn_Le        = 200.e+3        !  lateral diffusive length   [m]   (nn_aht_ijk_t= 0, 10) 
     
    890890   ln_dynvor_eeT = .false. !  energy conserving scheme (een using e3t) 
    891891   ln_dynvor_een = .false. !  energy & enstrophy scheme 
    892       nn_een_e3f = 0          ! =0  e3f = mi(mj(e3t))/4  
     892      nn_een_e3f = 0          ! =0  e3f = mi(mj(e3t))/4 
    893893      !                       ! =1  e3f = mi(mj(e3t))/mi(mj( tmask)) 
    894894   ln_dynvor_msk = .false. !  vorticity multiplied by fmask (=T)        ==>>> PLEASE DO NOT ACTIVATE 
     
    935935      !                             !  =-30  read in eddy_viscosity_3D.nc file 
    936936      !                             !  =-20  read in eddy_viscosity_2D.nc file 
    937       !                             !  =  0  constant  
     937      !                             !  =  0  constant 
    938938      !                             !  = 10  F(k)=c1d 
    939939      !                             !  = 20  F(i,j)=F(grid spacing)=c2d 
     
    941941      !                             !  = 31  F(i,j,k)=F(grid spacing and local velocity) 
    942942      !                             !  = 32  F(i,j,k)=F(local gridscale and deformation rate) 
    943       !                        !  time invariant coefficients :  ahm = 1/2  Uv*Lv   (lap case)  
     943      !                        !  time invariant coefficients :  ahm = 1/2  Uv*Lv   (lap case) 
    944944      !                             !                            or  = 1/12 Uv*Lv^3 (blp case) 
    945945      rn_Uv      = 0.1              !  lateral viscous velocity [m/s] (nn_ahm_ijk_t= 0, 10, 20, 30) 
     
    10651065                              !        = 0  constant 10 m length scale 
    10661066                              !        = 1  0.5m at the equator to 30m poleward of 40 degrees 
    1067       rn_eice     =   4       !  below sea ice: =0 ON ; =4 OFF when ice fraction > 1/4    
     1067      rn_eice     =   4       !  below sea ice: =0 ON ; =4 OFF when ice fraction > 1/4 
    10681068/ 
    10691069!----------------------------------------------------------------------- 
     
    11161116   rn_osm_mle_rho_c =    0.01      ! delta rho criterion used to calculate MLD for FK 
    11171117   rn_osm_mle_thresh = 0.0005      ! delta b criterion used for FK MLD criterion 
    1118    rn_osm_mle_tau     = 172800.   ! time scale for FK-OSM (seconds)  (case rn_osm_mle=0) 
    1119 / 
     1118   rn_osm_mle_tau     = 172800.    ! time scale for FK-OSM (seconds)  (case rn_osm_mle=0) 
     1119   ln_osm_hmle_limit   = .false.   ! limit hmle to rn_osm_hmle_limit*hbl 
     1120   rn_osm_hmle_limit   = 1.2 
     1121   / 
    11201122!----------------------------------------------------------------------- 
    11211123&namzdf_iwm    !    internal wave-driven mixing parameterization        (ln_zdfiwm =T) 
  • NEMO/branches/2019/dev_r11078_OSMOSIS_IMMERSE_Nurser_4.0/src/OCE/ZDF/zdfosm.F90

    r13396 r13398  
    123123   !                                        ! parameters used in nn_osm_mle = 1 case 
    124124   REAL(wp) ::   rn_osm_mle_lat              ! reference latitude for a 5 km scale of ML front 
     125   LOGICAL  ::   ln_osm_hmle_limit           ! If true arbitrarily restrict hmle to rn_osm_hmle_limit*zmld 
     126   REAL(wp) ::   rn_osm_hmle_limit           ! If ln_osm_hmle_limit true arbitrarily restrict hmle to rn_osm_hmle_limit*zmld 
    125127   REAL(wp) ::   rn_osm_mle_rho_c        ! Density criterion for definition of MLD used by FK 
    126128   REAL(wp) ::   r5_21 = 5.e0 / 21.e0   ! factor used in mle streamfunction computation 
     
    10701072 
    10711073 
    1072   
     1074 
    10731075       IF ( ln_osm_mle ) THEN  ! set up diffusivity and non-gradient mixing 
    10741076          DO jj = 2 , jpjm1 
     
    18211823      zN2_c = grav * rn_osm_mle_rho_c * r1_rau0   ! convert density criteria into N^2 criteria 
    18221824      DO jk = nlb10, jpkm1 
    1823          DO jj = 1, jpj                ! Mixed layer level: w-level  
     1825         DO jj = 1, jpj                ! Mixed layer level: w-level 
    18241826            DO ji = 1, jpi 
    18251827               ikt = mbkt(ji,jj) 
     
    19401942          ENDIF 
    19411943          hmle(ji,jj) = MIN(hmle(ji,jj), ht_n(ji,jj)) 
     1944          IF(ln_osm_hmle_limit) hmle(ji,jj) = MIN(hmle(ji,jj), rn_osm_hmle_limit*zmld(ji,jj)) 
    19421945        END DO 
    19431946      END DO 
     
    19701973             zwb_fk(ji,jj) = rn_osm_mle_ce * hmle(ji,jj) * hmle(ji,jj) *ztmp * zdbds_mle * zdbds_mle 
    19711974      ! This velocity scale, defined in Fox-Kemper et al (2008), is needed for calculating dhdt. 
    1972              zvel_mle(ji,jj) = zdbds_mle * ztmp * hmle(ji,jj) * tmask(ji,jj,1)  
     1975             zvel_mle(ji,jj) = zdbds_mle * ztmp * hmle(ji,jj) * tmask(ji,jj,1) 
    19731976             zdiff_mle(ji,jj) = 1.e-4_wp * ztmp * zdbds_mle * zhmle(ji,jj)**3  / rn_osm_mle_lf 
    19741977          ENDIF 
     
    20022005! Namelist for Fox-Kemper parametrization. 
    20032006      NAMELIST/namosm_mle/ nn_osm_mle, rn_osm_mle_ce, rn_osm_mle_lf, rn_osm_mle_time, rn_osm_mle_lat,& 
    2004            & rn_osm_mle_rho_c,rn_osm_mle_thresh,rn_osm_mle_tau 
     2007           & rn_osm_mle_rho_c, rn_osm_mle_thresh, rn_osm_mle_tau, ln_osm_hmle_limit, rn_osm_hmle_limit 
    20052008 
    20062009     !!---------------------------------------------------------------------- 
     
    20712074            WRITE(numout,*) '         Threshold used to define ML for FK                      rn_osm_mle_thresh  = ', rn_osm_mle_thresh, 'm^2/s' 
    20722075            WRITE(numout,*) '         Timescale for OSM-FK                                         rn_osm_mle_tau  = ', rn_osm_mle_tau, 's' 
     2076            WRITE(numout,*) '         switch to limit hmle                                      ln_osm_hmle_limit  = ', ln_osm_hmle_limit 
     2077            WRITE(numout,*) '         fraction of zmld to limit hmle to if ln_osm_hmle_limit =.T.  rn_osm_hmle_limit = ', rn_osm_hmle_limit 
    20732078         ENDIF         ! 
    20742079     ENDIF 
Note: See TracChangeset for help on using the changeset viewer.