Changeset 13883
- Timestamp:
- 2020-11-26T11:27:11+01:00 (4 years ago)
- Location:
- NEMO/branches/UKMO/NEMO_4.0.3_GM_rossby_radius_ramp
- Files:
-
- 4 edited
Legend:
- Unmodified
- Added
- Removed
-
NEMO/branches/UKMO/NEMO_4.0.3_GM_rossby_radius_ramp/cfgs/SHARED/field_def_nemo-oce.xml
r13340 r13883 92 92 <field id="beta" long_name="haline contraction" unit="1e3" grid_ref="grid_T_3D" /> 93 93 <field id="rhop" long_name="potential density (sigma0)" standard_name="sea_water_sigma_theta" unit="kg/m3" grid_ref="grid_T_3D" /> 94 <field id="N_2d" long_name="Depth-mean N" unit="m/s" /> 95 <field id="modslp" long_name="sqrt( slpi^2 + slpj^2 )" unit="1" grid_ref="grid_T_3D" /> 96 <field id="RossRad" long_name="Rossby radius" unit="m" /> 97 <field id="RossRadlim" long_name="Rossby radius (limited)" unit="m" /> 98 <field id="Tclinic_recip" long_name="recip of baroclinic timescale" unit="s-1" /> 99 <field id="RR_GS" long_name="Rossby radius / min(dx,dy)" unit="1" /> 94 100 95 101 <!-- Energy - horizontal divergence --> -
NEMO/branches/UKMO/NEMO_4.0.3_GM_rossby_radius_ramp/cfgs/SHARED/namelist_ref
r13587 r13883 837 837 ! ! = 30 F(i,j,k) =ldf_c2d * ldf_c1d 838 838 ! ! time invariant coefficients: aei0 = 1/2 Ue*Le 839 rn_Ue = 0.02 ! lateral diffusive velocity [m/s] (nn_a ht_ijk_t= 0, 10, 20, 30)840 rn_Le = 200.e+3 ! lateral diffusive length [m] (nn_a ht_ijk_t= 0, 10)839 rn_Ue = 0.02 ! lateral diffusive velocity [m/s] (nn_aei_ijk_t= 0, 10, 20, 30) 840 rn_Le = 200.e+3 ! lateral diffusive length [m] (nn_aei_ijk_t= 0, 10) 841 841 ! 842 842 ln_ldfeiv_dia =.false. ! diagnose eiv stream function and velocities 843 nn_ldfeiv_shape = 0 ! shape of bounding coefficient (nn_aei_ijk_t= 21 only) 843 844 / 844 845 !----------------------------------------------------------------------- -
NEMO/branches/UKMO/NEMO_4.0.3_GM_rossby_radius_ramp/src/OCE/DIA/diawri.F90
r13587 r13883 391 391 392 392 CALL iom_put( "bn2", rn2 ) ! Brunt-Vaisala buoyancy frequency (N^2) 393 IF( (.NOT.l_ldfeiv_time) .AND. ( iom_use('RossRad') .OR. iom_use('RossRadlim') & 394 & .OR. iom_use('TclinicR') .OR. iom_use('RR_GG') & 395 & .OR. iom_use('aeiu_2d') .OR. iom_use('aeiv_2d') ) ) THEN 396 CALL ldf_eiv(kt, 75.0, z2d, z3d(:,:,1)) 397 CALL iom_put('aeiu_2d', z2d) 398 CALL iom_put('aeiv_2d', z3d(:,:,1)) 399 ENDIF 393 400 ! 394 401 -
NEMO/branches/UKMO/NEMO_4.0.3_GM_rossby_radius_ramp/src/OCE/LDF/ldftra.F90
r13587 r13883 72 72 REAL(wp), PUBLIC :: rn_Ue !: lateral diffusive velocity [m/s] 73 73 REAL(wp), PUBLIC :: rn_Le !: lateral diffusive length [m] 74 INTEGER, PUBLIC :: nn_ldfeiv_shape !: shape of bounding coefficient (Treguier et al formulation only) 74 75 75 76 ! ! Flag to control the type of lateral diffusive operator … … 409 410 !!---------------------------------------------------------------------- 410 411 ! 411 IF( ln_ldfeiv .AND. nn_aei_ijk_t == 21 ) THEN ! eddy induced velocity coefficients 412 IF( ln_ldfeiv .AND. ( nn_aei_ijk_t == 21 ) ) THEN 413 ! ! eddy induced velocity coefficients 412 414 ! ! =F(growth rate of baroclinic instability) 413 415 ! ! max value aeiv_0 ; decreased to 0 within 20N-20S … … 502 504 !! 503 505 NAMELIST/namtra_eiv/ ln_ldfeiv , ln_ldfeiv_dia, & ! eddy induced velocity (eiv) 504 & nn_aei_ijk_t, rn_Ue, rn_Le ! eiv coefficient 506 & nn_aei_ijk_t, rn_Ue, rn_Le, & ! eiv coefficient 507 & nn_ldfeiv_shape 505 508 !!---------------------------------------------------------------------- 506 509 ! … … 525 528 WRITE(numout,*) ' eiv streamfunction & velocity diag. ln_ldfeiv_dia = ', ln_ldfeiv_dia 526 529 WRITE(numout,*) ' coefficients :' 527 WRITE(numout,*) ' type of time-space variation nn_aei_ijk_t = ', nn_a ht_ijk_t530 WRITE(numout,*) ' type of time-space variation nn_aei_ijk_t = ', nn_aei_ijk_t 528 531 WRITE(numout,*) ' lateral diffusive velocity (if cst) rn_Ue = ', rn_Ue, ' m/s' 529 532 WRITE(numout,*) ' lateral diffusive length (if cst) rn_Le = ', rn_Le, ' m' … … 595 598 IF(lwp) WRITE(numout,*) ' = F( growth rate of baroclinic instability )' 596 599 IF(lwp) WRITE(numout,*) ' maximum allowed value: aei0 = ', aei0, ' m2/s' 600 IF(lwp) WRITE(numout,*) ' shape of bounding coefficient : ',nn_ldfeiv_shape 597 601 ! 598 602 l_ldfeiv_time = .TRUE. ! will be calculated by call to ldf_tra routine in step.F90 … … 638 642 !!---------------------------------------------------------------------- 639 643 INTEGER , INTENT(in ) :: kt ! ocean time-step index 640 REAL(wp) , INTENT(in out) :: paei0 ! max value [m2/s]644 REAL(wp) , INTENT(in ) :: paei0 ! max value [m2/s] 641 645 REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) :: paeiu, paeiv ! eiv coefficient [m2/s] 642 646 ! 643 647 INTEGER :: ji, jj, jk ! dummy loop indices 644 REAL(wp) :: zfw, ze3w, zn2, z1_f20, zaht, zaht_min, zzaei ! local scalars 645 REAL(wp), DIMENSION(jpi,jpj) :: zn, zah, zhw, zRo, zaeiw ! 2D workspace 648 REAL(wp) :: zfw, ze3w, zn2, z1_f20, zaht, zaht_min, zzaei, z2_3 ! local scalars 649 REAL(wp), DIMENSION(jpi,jpj) :: zn, zah, zhw, zRo, zRo_lim, zTclinic_recip, zaeiw, zratio ! 2D workspace 650 REAL(wp), DIMENSION(jpi,jpj,jpk) :: zmodslp ! 3D workspace 646 651 !!---------------------------------------------------------------------- 647 652 ! 648 653 zn (:,:) = 0._wp ! Local initialization 654 zmodslp(:,:,:) = 0._wp 649 655 zhw(:,:) = 5._wp 650 656 zah(:,:) = 0._wp 651 657 zRo(:,:) = 0._wp 658 zRo_lim(:,:) = 0._wp 659 zTclinic_recip(:,:) = 0._wp 660 zratio(:,:) = 0._wp 661 zaeiw(:,:) = 0._wp 652 662 ! ! Compute lateral diffusive coefficient at T-point 653 663 IF( ln_traldf_triad ) THEN … … 682 692 ! T^-1 = sqrt(m_jpk(N^2*(r1^2+r2^2)*e3w)) 683 693 ze3w = e3w_n(ji,jj,jk) * wmask(ji,jj,jk) 684 zah(ji,jj) = zah(ji,jj) + zn2 * ( wslpi(ji,jj,jk) * wslpi(ji,jj,jk) & 685 & + wslpj(ji,jj,jk) * wslpj(ji,jj,jk) ) * ze3w 694 zmodslp(ji,jj,jk) = wslpi(ji,jj,jk) * wslpi(ji,jj,jk) & 695 & + wslpj(ji,jj,jk) * wslpj(ji,jj,jk) 696 zah(ji,jj) = zah(ji,jj) + zn2 * zmodslp(ji,jj,jk) * ze3w 686 697 zhw(ji,jj) = zhw(ji,jj) + ze3w 687 698 END DO … … 693 704 DO ji = fs_2, fs_jpim1 ! vector opt. 694 705 zfw = MAX( ABS( 2. * omega * SIN( rad * gphit(ji,jj) ) ) , 1.e-10 ) 695 ! Rossby radius at w-point taken betwenn 2 km and 40km 696 zRo(ji,jj) = MAX( 2.e3 , MIN( .4 * zn(ji,jj) / zfw, 40.e3 ) ) 706 ! Rossby radius at w-point taken between 2 km and 40km 707 zRo(ji,jj) = .4 * zn(ji,jj) / zfw 708 zRo_lim(ji,jj) = MAX( 2.e3 , MIN( zRo(ji,jj), 40.e3 ) ) 697 709 ! Compute aeiw by multiplying Ro^2 and T^-1 698 zaeiw(ji,jj) = zRo(ji,jj) * zRo(ji,jj) * SQRT( zah(ji,jj) / zhw(ji,jj) ) * tmask(ji,jj,1) 710 zTclinic_recip(ji,jj) = SQRT( MAX(zah(ji,jj),0._wp) / zhw(ji,jj) ) * tmask(ji,jj,1) 711 zaeiw(ji,jj) = zRo_lim(ji,jj) * zRo_lim(ji,jj) * zTclinic_recip(ji,jj) 699 712 END DO 700 713 END DO 701 714 IF( iom_use('N_2d') ) CALL iom_put('N_2d',zn(:,:)/zhw(:,:)) 715 IF( iom_use('modslp') ) CALL iom_put('modslp',SQRT(zmodslp(:,:,:)) ) 716 CALL iom_put('RossRad',zRo) 717 CALL iom_put('RossRadlim',zRo_lim) 718 CALL iom_put('Tclinic_recip',zTclinic_recip) 702 719 ! !== Bound on eiv coeff. ==! 703 720 z1_f20 = 1._wp / ( 2._wp * omega * sin( rad * 20._wp ) ) 704 DO jj = 2, jpjm1 705 DO ji = fs_2, fs_jpim1 ! vector opt. 706 zzaei = MIN( 1._wp, ABS( ff_t(ji,jj) * z1_f20 ) ) * zaeiw(ji,jj) ! tropical decrease 707 zaeiw(ji,jj) = MIN( zzaei , paei0 ) ! Max value = paei0 708 END DO 709 END DO 721 z2_3 = 2._wp/3._wp 722 723 SELECT CASE(nn_ldfeiv_shape) 724 CASE(0) !! Standard shape applied - decrease in tropics and cap. 725 DO jj = 2, jpjm1 726 DO ji = fs_2, fs_jpim1 ! vector opt. 727 zzaei = MIN( 1._wp, ABS( ff_t(ji,jj) * z1_f20 ) ) * zaeiw(ji,jj) ! tropical decrease 728 zaeiw(ji,jj) = MIN( zzaei, paei0 ) 729 END DO 730 END DO 731 732 CASE(1) !! Abrupt cut-off on Rossby radius: 733 ! JD : modifications here to introduce scaling by local rossby radius of deformation vs local grid scale 734 ! arbitrary decision that GM is de-activated if local rossy radius larger than 2 times local grid scale 735 ! based on Hallberg (2013) 736 DO jj = 2, jpjm1 737 DO ji = fs_2, fs_jpim1 ! vector opt. 738 IF ( zRo(ji,jj) >= ( 2._wp * MIN( e1t(ji,jj), e2t(ji,jj) ) ) ) THEN 739 ! TODO : use a version of zRo that integrates over a few time steps ? 740 zaeiw(ji,jj) = 0._wp 741 ELSE 742 zaeiw(ji,jj) = MIN( zaeiw(ji,jj), paei0 ) 743 ENDIF 744 END DO 745 END DO 746 747 CASE(2) !! Rossby radius ramp type 1: 748 DO jj = 2, jpjm1 749 DO ji = fs_2, fs_jpim1 ! vector opt. 750 zratio(ji,jj) = zRo(ji,jj)/MIN(e1t(ji,jj),e2t(ji,jj)) 751 zaeiw(ji,jj) = MIN( zaeiw(ji,jj), MAX( 0._wp, MIN( 1._wp, z2_3*(2._wp - zratio(ji,jj)) ) ) * paei0 ) 752 END DO 753 END DO 754 CALL iom_put('RR_GS',zratio) 755 756 CASE(3) !! Rossby radius ramp type 2: 757 DO jj = 2, jpjm1 758 DO ji = fs_2, fs_jpim1 ! vector opt. 759 zratio(ji,jj) = MIN(e1t(ji,jj),e2t(ji,jj))/zRo(ji,jj) 760 zaeiw(ji,jj) = MIN( zaeiw(ji,jj), MAX( 0._wp, MIN( 1._wp, z2_3*( zratio(ji,jj) - 0.5_wp ) ) ) * paei0 ) 761 END DO 762 END DO 763 764 CASE(4) !! bathymetry ramp: 765 DO jj = 2, jpjm1 766 DO ji = fs_2, fs_jpim1 ! vector opt. 767 zaeiw(ji,jj) = MIN( zaeiw(ji,jj), MAX( 0._wp, MIN( 1._wp, 0.001*(ht_0(ji,jj) - 2000._wp) ) ) * paei0 ) 768 END DO 769 END DO 770 771 CASE(5) !! Rossby radius ramp type 1 applied to Treguier et al coefficient rather than cap: 772 !! Note the ramp is RR/GS=[2.0,1.0] (not [2.0,0.5] as for cases 2,3) and we ramp up 773 !! to 5% of the Treguier et al coefficient, aiming for peak values of around 100m2/s 774 !! at high latitudes rather than 2000m2/s which is what you get in eORCA025 with an 775 !! uncapped coefficient. 776 DO jj = 2, jpjm1 777 DO ji = fs_2, fs_jpim1 ! vector opt. 778 zratio(ji,jj) = zRo(ji,jj)/MIN(e1t(ji,jj),e2t(ji,jj)) 779 zaeiw(ji,jj) = MAX( 0._wp, MIN( 1._wp, 2._wp - zratio(ji,jj) ) ) * 0.05 * zaeiw(ji,jj) 780 zaeiw(ji,jj) = MIN( zaeiw(ji,jj), paei0 ) 781 END DO 782 END DO 783 CALL iom_put('RR_GS',zratio) 784 785 CASE DEFAULT 786 CALL ctl_stop('ldf_eiv: Unrecognised option for nn_ldfeiv_shape.') 787 788 END SELECT 789 710 790 CALL lbc_lnk( 'ldftra', zaeiw(:,:), 'W', 1. ) ! lateral boundary condition 711 791 ! … … 722 802 paeiv(:,:,jk) = paeiv(:,:,1) * vmask(:,:,jk) 723 803 END DO 724 ! 804 ! 725 805 END SUBROUTINE ldf_eiv 726 806
Note: See TracChangeset
for help on using the changeset viewer.