[3] | 1 | MODULE zdfmxl |
---|
| 2 | !!====================================================================== |
---|
| 3 | !! *** MODULE zdfmxl *** |
---|
| 4 | !! Ocean physics: mixed layer depth |
---|
| 5 | !!====================================================================== |
---|
[463] | 6 | !! History : |
---|
| 7 | !! 9.0 ! 03-08 (G. Madec) OpenMP/autotasking optimization |
---|
[3] | 8 | !!---------------------------------------------------------------------- |
---|
| 9 | !! zdf_mxl : Compute the turbocline and mixed layer depths. |
---|
| 10 | !!---------------------------------------------------------------------- |
---|
| 11 | !! * Modules used |
---|
| 12 | USE oce ! ocean dynamics and tracers variables |
---|
| 13 | USE dom_oce ! ocean space and time domain variables |
---|
| 14 | USE zdf_oce ! ocean vertical physics |
---|
| 15 | USE in_out_manager ! I/O manager |
---|
[258] | 16 | USE prtctl ! Print control |
---|
[3] | 17 | |
---|
| 18 | IMPLICIT NONE |
---|
| 19 | PRIVATE |
---|
| 20 | |
---|
| 21 | !! * Routine accessibility |
---|
| 22 | PUBLIC zdf_mxl ! called by step.F90 |
---|
| 23 | |
---|
| 24 | !! * Shared module variables |
---|
[16] | 25 | INTEGER, PUBLIC, DIMENSION(jpi,jpj) :: & !: |
---|
[3] | 26 | nmln !: number of level in the mixed layer |
---|
[16] | 27 | REAL(wp), PUBLIC, DIMENSION(jpi,jpj) :: & !: |
---|
[3] | 28 | hmld , & !: mixing layer depth (turbocline) (m) |
---|
| 29 | hmlp , & !: mixed layer depth (rho=rho0+zdcrit) (m) |
---|
| 30 | hmlpt !: mixed layer depth at t-points (m) |
---|
| 31 | |
---|
| 32 | !! * module variables |
---|
| 33 | REAL(wp) :: & |
---|
| 34 | avt_c = 5.e-4_wp, & ! Kz criterion for the turbocline depth |
---|
| 35 | rho_c = 0.01_wp ! density criterion for mixed layer depth |
---|
| 36 | |
---|
| 37 | !! * Substitutions |
---|
| 38 | # include "domzgr_substitute.h90" |
---|
| 39 | !!---------------------------------------------------------------------- |
---|
[247] | 40 | !! OPA 9.0 , LOCEAN-IPSL (2005) |
---|
| 41 | !! $Header$ |
---|
| 42 | !! This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt |
---|
[3] | 43 | !!---------------------------------------------------------------------- |
---|
| 44 | |
---|
| 45 | CONTAINS |
---|
| 46 | |
---|
[463] | 47 | # if defined key_mpp_omp |
---|
[3] | 48 | !!---------------------------------------------------------------------- |
---|
[463] | 49 | !! 'key_mpp_omp' j-k-i loop (j-slab) |
---|
[3] | 50 | !!---------------------------------------------------------------------- |
---|
| 51 | |
---|
| 52 | SUBROUTINE zdf_mxl( kt ) |
---|
| 53 | !!---------------------------------------------------------------------- |
---|
| 54 | !! *** ROUTINE zdfmxl *** |
---|
| 55 | !! |
---|
| 56 | !! ** Purpose : Compute the turbocline depth and the mixed layer depth |
---|
| 57 | !! with a density criteria. |
---|
| 58 | !! |
---|
| 59 | !! ** Method : The turbocline depth is the depth at which the vertical |
---|
| 60 | !! eddy diffusivity coefficient (resulting from the vertical physics |
---|
| 61 | !! alone, not the isopycnal part, see trazdf.F) fall below a given |
---|
| 62 | !! value defined locally (avt_c here taken equal to 5 cm/s2) |
---|
| 63 | !! |
---|
| 64 | !! ** Action : |
---|
| 65 | !! |
---|
| 66 | !!---------------------------------------------------------------------- |
---|
| 67 | !! * Arguments |
---|
| 68 | INTEGER, INTENT( in ) :: kt ! ocean time-step index |
---|
| 69 | |
---|
| 70 | !! * Local declarations |
---|
| 71 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 72 | INTEGER :: ik ! temporary integer |
---|
| 73 | INTEGER, DIMENSION(jpi,jpj) :: & |
---|
| 74 | imld ! temporary workspace |
---|
| 75 | !!---------------------------------------------------------------------- |
---|
| 76 | |
---|
| 77 | IF( kt == nit000 ) THEN |
---|
| 78 | IF(lwp) WRITE(numout,*) |
---|
[463] | 79 | IF(lwp) WRITE(numout,*) 'zdf_mxl : mixed layer depth, j-k-i loops' |
---|
| 80 | IF(lwp) WRITE(numout,*) '~~~~~~~' |
---|
[3] | 81 | ENDIF |
---|
| 82 | |
---|
| 83 | ! ! =============== |
---|
| 84 | DO jj = 1, jpj ! Vertical slab |
---|
| 85 | ! ! =============== |
---|
| 86 | |
---|
| 87 | ! 1. Turbocline depth |
---|
| 88 | ! ------------------- |
---|
| 89 | ! last w-level at which avt<avt_c (starting from the bottom jk=jpk) |
---|
| 90 | ! (since avt(.,.,jpk)=0, we have jpk=< imld =< 2 ) |
---|
| 91 | DO jk = jpk, 2, -1 |
---|
| 92 | DO ji = 1, jpi |
---|
| 93 | IF( avt(ji,jj,jk) < avt_c ) imld(ji,jj) = jk |
---|
| 94 | END DO |
---|
| 95 | END DO |
---|
| 96 | |
---|
| 97 | ! Turbocline depth and sub-turbocline temperature |
---|
| 98 | DO ji = 1, jpi |
---|
| 99 | ik = imld(ji,jj) |
---|
| 100 | hmld (ji,jj) = fsdepw(ji,jj,ik) * tmask(ji,jj,1) |
---|
| 101 | END DO |
---|
| 102 | |
---|
| 103 | !!gm idea |
---|
| 104 | !! |
---|
| 105 | !!gm DO jk = jpk, 2, -1 |
---|
| 106 | !!gm DO ji = 1, jpi |
---|
| 107 | !!gm IF( avt(ji,jj,jk) < avt_c ) hmld(ji,jj) = fsdepw(ji,jj,jk) * tmask(ji,jj,1) |
---|
| 108 | !!gm END DO |
---|
| 109 | !!gm END DO |
---|
| 110 | !!gm |
---|
| 111 | |
---|
| 112 | ! 2. Mixed layer depth |
---|
| 113 | ! -------------------- |
---|
| 114 | ! Initialization to the number of w ocean point mbathy |
---|
| 115 | nmln(:,jj) = mbathy(:,jj) |
---|
| 116 | |
---|
| 117 | ! Last w-level at which rhop>=rho surf+rho_c (starting from jpk-1) |
---|
| 118 | ! (rhop defined at t-point, thus jk-1 for w-level just above) |
---|
| 119 | DO jk = jpkm1, 2, -1 |
---|
| 120 | DO ji = 1, jpi |
---|
| 121 | IF( rhop(ji,jj,jk) > rhop(ji,jj,1) + rho_c ) nmln(ji,jj) = jk |
---|
| 122 | END DO |
---|
| 123 | END DO |
---|
| 124 | |
---|
| 125 | ! Mixed layer depth |
---|
| 126 | DO ji = 1, jpi |
---|
| 127 | ik = nmln(ji,jj) |
---|
| 128 | hmlp (ji,jj) = fsdepw(ji,jj,ik) * tmask(ji,jj,1) |
---|
| 129 | hmlpt(ji,jj) = fsdept(ji,jj,ik-1) |
---|
| 130 | END DO |
---|
| 131 | ! ! =============== |
---|
| 132 | END DO ! End of slab |
---|
| 133 | ! ! =============== |
---|
| 134 | |
---|
[463] | 135 | IF(ln_ctl) CALL prt_ctl( tab2d_1=REAL(nmln,wp), clinfo1=' nmln : ', tab2d_2=hmld, clinfo2=' hmld : ', ovlap=1 ) |
---|
[49] | 136 | |
---|
[3] | 137 | END SUBROUTINE zdf_mxl |
---|
| 138 | |
---|
| 139 | # else |
---|
| 140 | !!---------------------------------------------------------------------- |
---|
| 141 | !! Default option : k-j-i loop |
---|
| 142 | !!---------------------------------------------------------------------- |
---|
| 143 | |
---|
| 144 | SUBROUTINE zdf_mxl( kt ) |
---|
| 145 | !!---------------------------------------------------------------------- |
---|
| 146 | !! *** ROUTINE zdfmxl *** |
---|
| 147 | !! |
---|
| 148 | !! ** Purpose : Compute the turbocline depth and the mixed layer depth |
---|
| 149 | !! with density criteria. |
---|
| 150 | !! |
---|
| 151 | !! ** Method : The turbocline depth is the depth at which the vertical |
---|
| 152 | !! eddy diffusivity coefficient (resulting from the vertical physics |
---|
| 153 | !! alone, not the isopycnal part, see trazdf.F) fall below a given |
---|
| 154 | !! value defined locally (avt_c here taken equal to 5 cm/s2) |
---|
| 155 | !! |
---|
| 156 | !! ** Action : |
---|
| 157 | !! |
---|
| 158 | !! History : |
---|
| 159 | !! ! 94-11 (M. Imbard) Original code |
---|
| 160 | !! 8.0 ! 96-01 (E. Guilyardi) sub mixed layer temp. |
---|
| 161 | !! 8.1 ! 97-07 (G. Madec) optimization |
---|
| 162 | !! 8.5 ! 02-06 (G. Madec) F90: Free form and module |
---|
| 163 | !!---------------------------------------------------------------------- |
---|
| 164 | !! * Arguments |
---|
| 165 | INTEGER, INTENT( in ) :: kt ! ocean time-step index |
---|
| 166 | |
---|
| 167 | !! * Local declarations |
---|
| 168 | INTEGER :: ji, jj, jk ! dummy loop indices |
---|
| 169 | INTEGER :: ik ! temporary integer |
---|
| 170 | INTEGER, DIMENSION(jpi,jpj) :: & |
---|
| 171 | imld ! temporary workspace |
---|
| 172 | !!---------------------------------------------------------------------- |
---|
| 173 | |
---|
| 174 | IF( kt == nit000 ) THEN |
---|
| 175 | IF(lwp) WRITE(numout,*) |
---|
| 176 | IF(lwp) WRITE(numout,*) 'zdf_mxl : mixed layer depth' |
---|
| 177 | IF(lwp) WRITE(numout,*) '~~~~~~~ ' |
---|
| 178 | ENDIF |
---|
| 179 | |
---|
| 180 | |
---|
| 181 | ! 1. Turbocline depth |
---|
| 182 | ! ------------------- |
---|
| 183 | ! last w-level at which avt<avt_c (starting from the bottom jk=jpk) |
---|
| 184 | ! (since avt(.,.,jpk)=0, we have jpk=< imld =< 2 ) |
---|
| 185 | DO jk = jpk, 2, -1 |
---|
| 186 | DO jj = 1, jpj |
---|
| 187 | DO ji = 1, jpi |
---|
| 188 | IF( avt(ji,jj,jk) < avt_c ) imld(ji,jj) = jk |
---|
| 189 | END DO |
---|
| 190 | END DO |
---|
| 191 | END DO |
---|
| 192 | |
---|
| 193 | ! Turbocline depth and sub-turbocline temperature |
---|
| 194 | DO jj = 1, jpj |
---|
| 195 | DO ji = 1, jpi |
---|
| 196 | ik = imld(ji,jj) |
---|
| 197 | hmld (ji,jj) = fsdepw(ji,jj,ik) * tmask(ji,jj,1) |
---|
| 198 | END DO |
---|
| 199 | END DO |
---|
| 200 | |
---|
| 201 | !!gm idea |
---|
| 202 | !! |
---|
| 203 | !!gm DO jk = jpk, 2, -1 |
---|
| 204 | !!gm DO jj = 1, jpj |
---|
| 205 | !!gm DO ji = 1, jpi |
---|
| 206 | !!gm IF( avt(ji,jj,jk) < avt_c ) hmld(ji,jj) = fsdepw(ji,jj,jk) * tmask(ji,jj,1) |
---|
| 207 | !!gm END DO |
---|
| 208 | !!gm END DO |
---|
| 209 | !!gm END DO |
---|
| 210 | !!gm |
---|
| 211 | |
---|
| 212 | ! 2. Mixed layer depth |
---|
| 213 | ! -------------------- |
---|
| 214 | ! Initialization to the number of w ocean point mbathy |
---|
| 215 | nmln(:,:) = mbathy(:,:) |
---|
| 216 | |
---|
| 217 | ! Last w-level at which rhop>=rho surf+rho_c (starting from jpk-1) |
---|
| 218 | ! (rhop defined at t-point, thus jk-1 for w-level just above) |
---|
| 219 | DO jk = jpkm1, 2, -1 |
---|
| 220 | DO jj = 1, jpj |
---|
| 221 | DO ji = 1, jpi |
---|
| 222 | IF( rhop(ji,jj,jk) > rhop(ji,jj,1) + rho_c ) nmln(ji,jj) = jk |
---|
| 223 | END DO |
---|
| 224 | END DO |
---|
| 225 | END DO |
---|
| 226 | |
---|
| 227 | ! Mixed layer depth |
---|
| 228 | DO jj = 1, jpj |
---|
| 229 | DO ji = 1, jpi |
---|
| 230 | ik = nmln(ji,jj) |
---|
| 231 | hmlp (ji,jj) = fsdepw(ji,jj,ik) * tmask(ji,jj,1) |
---|
| 232 | hmlpt(ji,jj) = fsdept(ji,jj,ik-1) |
---|
| 233 | END DO |
---|
| 234 | END DO |
---|
| 235 | |
---|
[463] | 236 | IF(ln_ctl) CALL prt_ctl( tab2d_1=REAL(nmln,wp), clinfo1=' nmln : ', tab2d_2=hmld, clinfo2=' hmld : ', ovlap=1 ) |
---|
[49] | 237 | |
---|
[3] | 238 | END SUBROUTINE zdf_mxl |
---|
| 239 | #endif |
---|
| 240 | |
---|
| 241 | !!====================================================================== |
---|
| 242 | END MODULE zdfmxl |
---|