Changes from NEMO/branches/2021/dev_r14393_HPC-03_Mele_Comm_Cleanup/src/OCE/ZDF/zdfmxl.F90 at r14757 to NEMO/branches/2021/dev_r14273_HPC-02_Daley_Tiling/src/OCE/ZDF/zdfmxl.F90 at r14780
- File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
NEMO/branches/2021/dev_r14273_HPC-02_Daley_Tiling/src/OCE/ZDF/zdfmxl.F90
r14757 r14780 26 26 PRIVATE 27 27 28 PUBLIC zdf_mxl ! called by zdfphy.F9028 PUBLIC zdf_mxl, zdf_mxl_turb ! called by zdfphy.F90 29 29 30 30 INTEGER , PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:) :: nmln !: number of level in the mixed layer (used by LDF, ZDF, TRD, TOP) … … 41 41 !!---------------------------------------------------------------------- 42 42 !! NEMO/OCE 4.0 , NEMO Consortium (2018) 43 !! $Id$ 43 !! $Id$ 44 44 !! Software governed by the CeCILL license (see ./LICENSE) 45 45 !!---------------------------------------------------------------------- … … 65 65 !! *** ROUTINE zdfmxl *** 66 66 !! 67 !! ** Purpose : Compute the turbocline depth and the mixed layer depth 68 !! with density criteria. 67 !! ** Purpose : Compute the mixed layer depth with density criteria. 69 68 !! 70 69 !! ** Method : The mixed layer depth is the shallowest W depth with 71 70 !! the density of the corresponding T point (just bellow) bellow a 72 71 !! given value defined locally as rho(10m) + rho_c 73 !! The turbocline depth is the depth at which the vertical74 !! eddy diffusivity coefficient (resulting from the vertical physics75 !! alone, not the isopycnal part, see trazdf.F) fall below a given76 !! value defined locally (avt_c here taken equal to 5 cm/s2 by default)77 72 !! 78 !! ** Action : nmln, hml d, hmlp, hmlpt73 !! ** Action : nmln, hmlp, hmlpt 79 74 !!---------------------------------------------------------------------- 80 75 INTEGER, INTENT(in) :: kt ! ocean time-step index … … 82 77 ! 83 78 INTEGER :: ji, jj, jk ! dummy loop indices 84 INTEGER :: iik n, iiki, ikt! local integer79 INTEGER :: iik, ikt ! local integer 85 80 REAL(wp) :: zN2_c ! local scalar 86 INTEGER, DIMENSION(jpi,jpj) :: imld ! 2D workspace87 81 !!---------------------------------------------------------------------- 88 82 ! 89 IF( kt == nit000 ) THEN 90 IF(lwp) WRITE(numout,*) 91 IF(lwp) WRITE(numout,*) 'zdf_mxl : mixed layer depth' 92 IF(lwp) WRITE(numout,*) '~~~~~~~ ' 93 ! ! allocate zdfmxl arrays 94 IF( zdf_mxl_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'zdf_mxl : unable to allocate arrays' ) 83 IF( .NOT. l_istiled .OR. ntile == 1 ) THEN ! Do only on the first tile 84 IF( kt == nit000 ) THEN 85 IF(lwp) WRITE(numout,*) 86 IF(lwp) WRITE(numout,*) 'zdf_mxl : mixed layer depth' 87 IF(lwp) WRITE(numout,*) '~~~~~~~ ' 88 ! ! allocate zdfmxl arrays 89 IF( zdf_mxl_alloc() /= 0 ) CALL ctl_stop( 'STOP', 'zdf_mxl : unable to allocate arrays' ) 90 ENDIF 95 91 ENDIF 96 92 ! 97 93 ! w-level of the mixing and mixed layers 98 nmln(:,:) = nlb10 ! Initialization to the number of w ocean point 99 hmlp(:,:) = 0._wp ! here hmlp used as a dummy variable, integrating vertically N^2 94 DO_2D_OVR( nn_hls, nn_hls, nn_hls, nn_hls ) 95 nmln(ji,jj) = nlb10 ! Initialization to the number of w ocean point 96 hmlp(ji,jj) = 0._wp ! here hmlp used as a dummy variable, integrating vertically N^2 97 END_2D 100 98 zN2_c = grav * rho_c * r1_rho0 ! convert density criteria into N^2 criteria 101 DO_3D ( nn_hls, nn_hls, nn_hls, nn_hls, nlb10, jpkm1 ) ! Mixed layer level: w-level99 DO_3D_OVR( nn_hls, nn_hls, nn_hls, nn_hls, nlb10, jpkm1 ) ! Mixed layer level: w-level 102 100 ikt = mbkt(ji,jj) 103 101 hmlp(ji,jj) = & … … 105 103 IF( hmlp(ji,jj) < zN2_c ) nmln(ji,jj) = MIN( jk , ikt ) + 1 ! Mixed layer level 106 104 END_3D 107 ! 108 ! w-level of the turbocline and mixing layer (iom_use) 109 imld(:,:) = mbkt(:,:) + 1 ! Initialization to the number of w ocean point 110 DO_3DS( nn_hls, nn_hls, nn_hls, nn_hls, jpkm1, nlb10, -1 ) ! from the bottom to nlb10 111 IF( avt (ji,jj,jk) < avt_c * wmask(ji,jj,jk) ) imld(ji,jj) = jk ! Turbocline 112 END_3D 113 ! depth of the mixing and mixed layers 114 DO_2D( nn_hls, nn_hls, nn_hls, nn_hls ) 115 iiki = imld(ji,jj) 116 iikn = nmln(ji,jj) 117 hmld (ji,jj) = gdepw(ji,jj,iiki ,Kmm) * ssmask(ji,jj) ! Turbocline depth 118 hmlp (ji,jj) = gdepw(ji,jj,iikn ,Kmm) * ssmask(ji,jj) ! Mixed layer depth 119 hmlpt(ji,jj) = gdept(ji,jj,iikn-1,Kmm) * ssmask(ji,jj) ! depth of the last T-point inside the mixed layer 105 ! depth of the mixed layer 106 DO_2D_OVR( nn_hls, nn_hls, nn_hls, nn_hls ) 107 iik = nmln(ji,jj) 108 hmlp (ji,jj) = gdepw(ji,jj,iik ,Kmm) * ssmask(ji,jj) ! Mixed layer depth 109 hmlpt(ji,jj) = gdept(ji,jj,iik-1,Kmm) * ssmask(ji,jj) ! depth of the last T-point inside the mixed layer 120 110 END_2D 121 111 ! 122 IF( .NOT.l_offline ) THEN 123 IF( iom_use("mldr10_1") ) THEN 124 IF( ln_isfcav ) THEN ; CALL iom_put( "mldr10_1", hmlp - risfdep) ! mixed layer thickness 125 ELSE ; CALL iom_put( "mldr10_1", hmlp ) ! mixed layer depth 126 END IF 112 IF( .NOT.l_offline .AND. iom_use("mldr10_1") ) THEN 113 IF( ln_isfcav ) THEN ; CALL iom_put( "mldr10_1", hmlp - risfdep) ! mixed layer thickness 114 ELSE ; CALL iom_put( "mldr10_1", hmlp ) ! mixed layer depth 127 115 END IF 128 IF( iom_use("mldkz5") ) THEN129 IF( ln_isfcav ) THEN ; CALL iom_put( "mldkz5" , hmld - risfdep ) ! turbocline thickness130 ELSE ; CALL iom_put( "mldkz5" , hmld ) ! turbocline depth131 END IF132 ENDIF133 116 ENDIF 134 117 ! … … 137 120 END SUBROUTINE zdf_mxl 138 121 122 123 SUBROUTINE zdf_mxl_turb( kt, Kmm ) 124 !!---------------------------------------------------------------------- 125 !! *** ROUTINE zdf_mxl_turb *** 126 !! 127 !! ** Purpose : Compute the turbocline depth. 128 !! 129 !! ** Method : The turbocline depth is the depth at which the vertical 130 !! eddy diffusivity coefficient (resulting from the vertical physics 131 !! alone, not the isopycnal part, see trazdf.F) fall below a given 132 !! value defined locally (avt_c here taken equal to 5 cm/s2 by default) 133 !! 134 !! ** Action : hmld 135 !!---------------------------------------------------------------------- 136 INTEGER, INTENT(in) :: kt ! ocean time-step index 137 INTEGER, INTENT(in) :: Kmm ! ocean time level index 138 ! 139 INTEGER :: ji, jj, jk ! dummy loop indices 140 INTEGER :: iik ! local integer 141 INTEGER, DIMENSION(A2D(nn_hls)) :: imld ! 2D workspace 142 !!---------------------------------------------------------------------- 143 ! 144 ! w-level of the turbocline and mixing layer (iom_use) 145 imld(:,:) = mbkt(A2D(nn_hls)) + 1 ! Initialization to the number of w ocean point 146 DO_3DS( 1, 1, 1, 1, jpkm1, nlb10, -1 ) ! from the bottom to nlb10 147 IF( avt (ji,jj,jk) < avt_c * wmask(ji,jj,jk) ) imld(ji,jj) = jk ! Turbocline 148 END_3D 149 ! depth of the mixing layer 150 DO_2D_OVR( 1, 1, 1, 1 ) 151 iik = imld(ji,jj) 152 hmld (ji,jj) = gdepw(ji,jj,iik ,Kmm) * ssmask(ji,jj) ! Turbocline depth 153 END_2D 154 ! 155 IF( .NOT.l_offline .AND. iom_use("mldkz5") ) THEN 156 IF( ln_isfcav ) THEN ; CALL iom_put( "mldkz5" , hmld - risfdep ) ! turbocline thickness 157 ELSE ; CALL iom_put( "mldkz5" , hmld ) ! turbocline depth 158 END IF 159 ENDIF 160 ! 161 END SUBROUTINE zdf_mxl_turb 139 162 !!====================================================================== 140 163 END MODULE zdfmxl
Note: See TracChangeset
for help on using the changeset viewer.