Changeset 7646 for trunk/NEMOGCM/NEMO/OPA_SRC/LDF/ldftra.F90
- Timestamp:
- 2017-02-06T10:25:03+01:00 (7 years ago)
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/NEMOGCM/NEMO/OPA_SRC/LDF/ldftra.F90
r6140 r7646 24 24 USE ldfslp ! lateral diffusion: slope of iso-neutral surfaces 25 25 USE ldfc1d_c2d ! lateral diffusion: 1D & 2D cases 26 USE dia ar5, ONLY: lk_diaar526 USE diaptr 27 27 ! 28 USE trc_oce, ONLY: lk_offline ! offline flag29 28 USE in_out_manager ! I/O manager 30 29 USE iom ! I/O module for ehanced bottom friction file … … 298 297 ! 299 298 INTEGER :: ji, jj, jk ! dummy loop indices 300 REAL(wp) :: zaht, zah t_min, z1_f20 ! local scalar301 !!---------------------------------------------------------------------- 302 ! 303 IF( nn_aei_ijk_t == 21 ) THEN ! eddy induced velocity coefficients299 REAL(wp) :: zaht, zahf, zaht_min, z1_f20 ! local scalar 300 !!---------------------------------------------------------------------- 301 ! 302 IF( ln_ldfeiv .AND. nn_aei_ijk_t == 21 ) THEN ! eddy induced velocity coefficients 304 303 ! ! =F(growth rate of baroclinic instability) 305 304 ! ! max value rn_aeiv_0 ; decreased to 0 within 20N-20S 306 305 CALL ldf_eiv( kt, rn_aeiv_0, aeiu, aeiv ) 307 IF(lwp .AND. kt<=nit000+20 ) WRITE(numout,*) ' kt , ldf_eiv appel', kt308 306 ENDIF 309 307 ! … … 314 312 ! ! max value rn_aht_0 (rn_aeiv_0 if nn_aei_ijk_t=21) 315 313 ! ! increase to rn_aht_0 within 20N-20S 316 IF( nn_aei_ijk_t /= 21 ) THEN 317 CALL ldf_eiv( kt, rn_aht_0, ahtu, ahtv ) 318 IF(lwp .AND. kt<=nit000+20 ) WRITE(numout,*) ' kt , ldf_eiv appel 2', kt 319 ELSE 314 IF( ln_ldfeiv .AND. nn_aei_ijk_t == 21 ) THEN ! use the already computed aei. 320 315 ahtu(:,:,1) = aeiu(:,:,1) 321 316 ahtv(:,:,1) = aeiv(:,:,1) 322 IF(lwp .AND. kt<=nit000+20 ) WRITE(numout,*) ' kt , ahtu=aeiu', kt 317 ELSE ! compute aht. 318 CALL ldf_eiv( kt, rn_aht_0, ahtu, ahtv ) 323 319 ENDIF 324 320 ! … … 327 323 DO jj = 1, jpj 328 324 DO ji = 1, jpi 329 zaht = ( 1._wp - MIN( 1._wp , ABS( ff(ji,jj) * z1_f20 ) ) ) * ( rn_aht_0 - zaht_min ) 325 !!gm CAUTION : here we assume lat/lon grid in 20deg N/S band (like all ORCA cfg) 326 !! ==>>> The Coriolis value is identical for t- & u_points, and for v- and f-points 327 zaht = ( 1._wp - MIN( 1._wp , ABS( ff_t(ji,jj) * z1_f20 ) ) ) * ( rn_aht_0 - zaht_min ) 328 zahf = ( 1._wp - MIN( 1._wp , ABS( ff_f(ji,jj) * z1_f20 ) ) ) * ( rn_aht_0 - zaht_min ) 330 329 ahtu(ji,jj,1) = ( MAX( zaht_min, ahtu(ji,jj,1) ) + zaht ) * umask(ji,jj,1) ! min value zaht_min 331 ahtv(ji,jj,1) = ( MAX( zaht_min, ahtv(ji,jj,1) ) + zah t) * vmask(ji,jj,1) ! increase within 20S-20N330 ahtv(ji,jj,1) = ( MAX( zaht_min, ahtv(ji,jj,1) ) + zahf ) * vmask(ji,jj,1) ! increase within 20S-20N 332 331 END DO 333 332 END DO … … 352 351 END SELECT 353 352 ! 354 IF( .NOT.lk_offline ) THEN 355 CALL iom_put( "ahtu_2d", ahtu(:,:,1) ) ! surface u-eddy diffusivity coeff. 356 CALL iom_put( "ahtv_2d", ahtv(:,:,1) ) ! surface v-eddy diffusivity coeff. 357 CALL iom_put( "ahtu_3d", ahtu(:,:,:) ) ! 3D u-eddy diffusivity coeff. 358 CALL iom_put( "ahtv_3d", ahtv(:,:,:) ) ! 3D v-eddy diffusivity coeff. 359 ! 353 CALL iom_put( "ahtu_2d", ahtu(:,:,1) ) ! surface u-eddy diffusivity coeff. 354 CALL iom_put( "ahtv_2d", ahtv(:,:,1) ) ! surface v-eddy diffusivity coeff. 355 CALL iom_put( "ahtu_3d", ahtu(:,:,:) ) ! 3D u-eddy diffusivity coeff. 356 CALL iom_put( "ahtv_3d", ahtv(:,:,:) ) ! 3D v-eddy diffusivity coeff. 357 ! 360 358 !!gm : THE IF below is to be checked (comes from Seb) 361 IF( ln_ldfeiv ) THEN 362 CALL iom_put( "aeiu_2d", aeiu(:,:,1) ) ! surface u-EIV coeff. 363 CALL iom_put( "aeiv_2d", aeiv(:,:,1) ) ! surface v-EIV coeff. 364 CALL iom_put( "aeiu_3d", aeiu(:,:,:) ) ! 3D u-EIV coeff. 365 CALL iom_put( "aeiv_3d", aeiv(:,:,:) ) ! 3D v-EIV coeff. 366 ENDIF 359 IF( ln_ldfeiv ) THEN 360 CALL iom_put( "aeiu_2d", aeiu(:,:,1) ) ! surface u-EIV coeff. 361 CALL iom_put( "aeiv_2d", aeiv(:,:,1) ) ! surface v-EIV coeff. 362 CALL iom_put( "aeiu_3d", aeiu(:,:,:) ) ! 3D u-EIV coeff. 363 CALL iom_put( "aeiv_3d", aeiv(:,:,:) ) ! 3D v-EIV coeff. 367 364 ENDIF 368 365 ! … … 555 552 END DO 556 553 557 !!gm IF( cp_cfg == "orca" .AND. jp_cfg == 2 ) THEN ! ORCA R2558 !!gm DO jj = 2, jpjm1559 !!gm DO ji = fs_2, fs_jpim1 ! vector opt.560 !!gm ! Take the minimum between aeiw and 1000 m2/s over shelves (depth shallower than 650 m)561 !!gm IF( mbkt(ji,jj) <= 20 ) zaeiw(ji,jj) = MIN( zaeiw(ji,jj), 1000. )562 !!gm END DO563 !!gm END DO564 !!gm ENDIF565 566 554 ! !== Bound on eiv coeff. ==! 567 555 z1_f20 = 1._wp / ( 2._wp * omega * sin( rad * 20._wp ) ) 568 556 DO jj = 2, jpjm1 569 557 DO ji = fs_2, fs_jpim1 ! vector opt. 570 zzaei = MIN( 1._wp, ABS( ff (ji,jj) * z1_f20 ) ) * zaeiw(ji,jj) ! tropical decrease558 zzaei = MIN( 1._wp, ABS( ff_t(ji,jj) * z1_f20 ) ) * zaeiw(ji,jj) ! tropical decrease 571 559 zaeiw(ji,jj) = MIN( zzaei , paei0 ) ! Max value = paei0 572 560 END DO … … 730 718 CALL iom_put( "woce_eiv", zw3d ) 731 719 ! 720 ! 721 ! 722 CALL wrk_alloc( jpi,jpj, zw2d ) 723 ! 724 zztmp = 0.5_wp * rau0 * rcp 725 IF( iom_use('ueiv_heattr') .OR. iom_use('ueiv_heattr3d') ) THEN 726 zw2d(:,:) = 0._wp 727 zw3d(:,:,:) = 0._wp 728 DO jk = 1, jpkm1 729 DO jj = 2, jpjm1 730 DO ji = fs_2, fs_jpim1 ! vector opt. 731 zw3d(ji,jj,jk) = zw3d(ji,jj,jk) + ( psi_uw(ji,jj,jk+1) - psi_uw(ji,jj,jk) ) & 732 & * ( tsn (ji,jj,jk,jp_tem) + tsn (ji+1,jj,jk,jp_tem) ) 733 zw2d(ji,jj) = zw2d(ji,jj) + zw3d(ji,jj,jk) 734 END DO 735 END DO 736 END DO 737 CALL lbc_lnk( zw2d, 'U', -1. ) 738 CALL lbc_lnk( zw3d, 'U', -1. ) 739 CALL iom_put( "ueiv_heattr" , zztmp * zw2d ) ! heat transport in i-direction 740 CALL iom_put( "ueiv_heattr3d", zztmp * zw3d ) ! heat transport in i-direction 741 ENDIF 742 zw2d(:,:) = 0._wp 743 zw3d(:,:,:) = 0._wp 744 DO jk = 1, jpkm1 745 DO jj = 2, jpjm1 746 DO ji = fs_2, fs_jpim1 ! vector opt. 747 zw3d(ji,jj,jk) = zw3d(ji,jj,jk) + ( psi_vw(ji,jj,jk+1) - psi_vw(ji,jj,jk) ) & 748 & * ( tsn (ji,jj,jk,jp_tem) + tsn (ji,jj+1,jk,jp_tem) ) 749 zw2d(ji,jj) = zw2d(ji,jj) + zw3d(ji,jj,jk) 750 END DO 751 END DO 752 END DO 753 CALL lbc_lnk( zw2d, 'V', -1. ) 754 CALL iom_put( "veiv_heattr", zztmp * zw2d ) ! heat transport in j-direction 755 CALL iom_put( "veiv_heattr", zztmp * zw3d ) ! heat transport in j-direction 756 ! 757 IF( ln_diaptr ) CALL dia_ptr_hst( jp_tem, 'eiv', 0.5 * zw3d ) 758 ! 759 zztmp = 0.5_wp * 0.5 760 IF( iom_use('ueiv_salttr') .OR. iom_use('ueiv_salttr3d')) THEN 761 zw2d(:,:) = 0._wp 762 zw3d(:,:,:) = 0._wp 763 DO jk = 1, jpkm1 764 DO jj = 2, jpjm1 765 DO ji = fs_2, fs_jpim1 ! vector opt. 766 zw3d(ji,jj,jk) = zw3d(ji,jj,jk) * ( psi_uw(ji,jj,jk+1) - psi_uw(ji,jj,jk) ) & 767 & * ( tsn (ji,jj,jk,jp_sal) + tsn (ji+1,jj,jk,jp_sal) ) 768 zw2d(ji,jj) = zw2d(ji,jj) + zw3d(ji,jj,jk) 769 END DO 770 END DO 771 END DO 772 CALL lbc_lnk( zw2d, 'U', -1. ) 773 CALL lbc_lnk( zw3d, 'U', -1. ) 774 CALL iom_put( "ueiv_salttr", zztmp * zw2d ) ! salt transport in i-direction 775 CALL iom_put( "ueiv_salttr3d", zztmp * zw3d ) ! salt transport in i-direction 776 ENDIF 777 zw2d(:,:) = 0._wp 778 zw3d(:,:,:) = 0._wp 779 DO jk = 1, jpkm1 780 DO jj = 2, jpjm1 781 DO ji = fs_2, fs_jpim1 ! vector opt. 782 zw3d(ji,jj,jk) = zw3d(ji,jj,jk) + ( psi_vw(ji,jj,jk+1) - psi_vw(ji,jj,jk) ) & 783 & * ( tsn (ji,jj,jk,jp_sal) + tsn (ji,jj+1,jk,jp_sal) ) 784 zw2d(ji,jj) = zw2d(ji,jj) + zw3d(ji,jj,jk) 785 END DO 786 END DO 787 END DO 788 CALL lbc_lnk( zw2d, 'V', -1. ) 789 CALL iom_put( "veiv_salttr", zztmp * zw2d ) ! salt transport in j-direction 790 CALL iom_put( "veiv_salttr", zztmp * zw3d ) ! salt transport in j-direction 791 ! 792 IF( ln_diaptr ) CALL dia_ptr_hst( jp_sal, 'eiv', 0.5 * zw3d ) 793 ! 794 CALL wrk_dealloc( jpi,jpj, zw2d ) 732 795 CALL wrk_dealloc( jpi,jpj,jpk, zw3d ) 733 !734 !735 IF( lk_diaar5 ) THEN !== eiv heat transport: calculate and output ==!736 CALL wrk_alloc( jpi,jpj, zw2d )737 !738 zztmp = 0.5_wp * rau0 * rcp739 zw2d(:,:) = 0._wp740 DO jk = 1, jpkm1741 DO jj = 2, jpjm1742 DO ji = fs_2, fs_jpim1 ! vector opt.743 zw2d(ji,jj) = zw2d(ji,jj) + zztmp * ( psi_uw(ji,jj,jk+1) - psi_uw(ji,jj,jk) ) &744 & * ( tsn (ji,jj,jk,jp_tem) + tsn (ji+1,jj,jk,jp_tem) )745 END DO746 END DO747 END DO748 CALL lbc_lnk( zw2d, 'U', -1. )749 CALL iom_put( "ueiv_heattr", zw2d ) ! heat transport in i-direction750 zw2d(:,:) = 0._wp751 DO jk = 1, jpkm1752 DO jj = 2, jpjm1753 DO ji = fs_2, fs_jpim1 ! vector opt.754 zw2d(ji,jj) = zw2d(ji,jj) + zztmp * ( psi_vw(ji,jj,jk+1) - psi_vw(ji,jj,jk) ) &755 & * ( tsn (ji,jj,jk,jp_tem) + tsn (ji,jj+1,jk,jp_tem) )756 END DO757 END DO758 END DO759 CALL lbc_lnk( zw2d, 'V', -1. )760 CALL iom_put( "veiv_heattr", zw2d ) ! heat transport in i-direction761 !762 CALL wrk_dealloc( jpi,jpj, zw2d )763 ENDIF764 796 ! 765 797 IF( nn_timing == 1 ) CALL timing_stop( 'ldf_eiv_dia')
Note: See TracChangeset
for help on using the changeset viewer.