/[MITgcm]/MITgcm/pkg/kpp/kpp_calc.F
ViewVC logotype

Contents of /MITgcm/pkg/kpp/kpp_calc.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.5 - (show annotations) (download)
Wed Nov 29 22:29:23 2000 UTC (23 years, 6 months ago) by adcroft
Branch: MAIN
CVS Tags: branch-atmos-merge-shapiro, branch-atmos-merge-freeze, branch-atmos-merge-start, checkpoint33, checkpoint34, branch-atmos-merge-zonalfilt, branch-atmos-merge-phase5, branch-atmos-merge-phase4, branch-atmos-merge-phase7, branch-atmos-merge-phase6, branch-atmos-merge-phase1, branch-atmos-merge-phase3, branch-atmos-merge-phase2
Branch point for: branch-atmos-merge
Changes since 1.4: +3 -2 lines
Fixed confusion about units of forcing arrays in FFIELDS.h
namely Fu,Fv,Qnet,Qsw,EmPmR:
  - Removed verification/natl_box/code/external_fields_scale.F
        (did not differ from that in model/src)
  - Changed units of fu,fv,Qnet,Qsw,EmPmR back to proper units
     (see FFIELDS.h for description)
  - Scale fu,fv,Qnet,Qsw,EmPmR when used in external_forcing_surf.F,
    kpp_calc.F and kpp_transport_t.F
  - Removed model/src/external_fields_scale.F and calls to it
  - verification/natl_box uses flux data with "atmospheric" sign so
    a special version of external_fields_load.F is used to
    change the data as it's read in. This way, the arrays
    have the right units and signs at all times tha a user could
    possibly use them.

1 C $Header: /u/gcmpack/models/MITgcmUV/pkg/kpp/kpp_calc.F,v 1.4 2000/11/13 16:37:02 heimbach Exp $
2
3 #include "KPP_OPTIONS.h"
4
5 subroutine KPP_CALC(
6 I bi, bj, myTime, myThid )
7 C /==========================================================\
8 C | SUBROUTINE KPP_CALC |
9 C | o Compute all KPP fields defined in KPP.h |
10 C |==========================================================|
11 C | This subroutine serves as an interface between MITGCMUV |
12 C | code and NCOM 1-D routines in kpp_routines.F |
13 C \==========================================================/
14 IMPLICIT NONE
15
16 c=======================================================================
17 c
18 c written by : jan morzel, august 11, 1994
19 c modified by : jan morzel, january 25, 1995 : "dVsq" and 1d code
20 c detlef stammer, august, 1997 : for MIT GCM Classic
21 c d. menemenlis, july, 1998 : for MIT GCM UV
22 c
23 c compute vertical mixing coefficients based on the k-profile
24 c and oceanic planetary boundary layer scheme by large & mcwilliams.
25 c
26 c summary:
27 c - compute interior mixing everywhere:
28 c interior mixing gets computed at all interfaces due to constant
29 c internal wave background activity ("fkpm" and "fkph"), which
30 c is enhanced in places of static instability (local richardson
31 c number < 0).
32 c Additionally, mixing can be enhanced by adding contribution due
33 c to shear instability which is a function of the local richardson
34 c number
35 c - double diffusivity:
36 c interior mixing can be enhanced by double diffusion due to salt
37 c fingering and diffusive convection (ifdef "kmixdd").
38 c - kpp scheme in the boundary layer:
39 c
40 c a.boundary layer depth:
41 c at every gridpoint the depth of the oceanic boundary layer
42 c ("hbl") gets computed by evaluating bulk richardson numbers.
43 c b.boundary layer mixing:
44 c within the boundary layer, above hbl, vertical mixing is
45 c determined by turbulent surface fluxes, and interior mixing at
46 c the lower boundary, i.e. at hbl.
47 c
48 c this subroutine provides the interface between the MIT GCM UV and the
49 c subroutine "kppmix", where boundary layer depth, vertical
50 c viscosity, vertical diffusivity, and counter gradient term (ghat)
51 c are computed slabwise.
52 c note: subroutine "kppmix" uses m-k-s units.
53 c
54 c time level:
55 c input tracer and velocity profiles are evaluated at time level
56 c tau, surface fluxes come from tau or tau-1.
57 c
58 c grid option:
59 c in this "1-grid" implementation, diffusivity and viscosity
60 c profiles are computed on the "t-grid" (by using velocity shear
61 c profiles averaged from the "u,v-grid" onto the "t-grid"; note, that
62 c the averaging includes zero values on coastal and seafloor grid
63 c points). viscosity on the "u,v-grid" is computed by averaging the
64 c "t-grid" viscosity values onto the "u,v-grid".
65 c
66 c vertical grid:
67 c mixing coefficients get evaluated at the bottom of the lowest
68 c layer, i.e., at depth zw(Nr). these values are only useful when
69 c the model ocean domain does not include the entire ocean down to
70 c the seafloor ("upperocean" setup) and allows flux through the
71 c bottom of the domain. for full-depth runs, these mixing
72 c coefficients are being zeroed out before leaving this subroutine.
73 c
74 c-------------------------------------------------------------------------
75
76 c global parameters updated by kpp_calc
77 c KPPviscAz - KPP eddy viscosity coefficient (m^2/s)
78 c KPPdiffKzT - KPP diffusion coefficient for temperature (m^2/s)
79 c KPPdiffKzS - KPP diffusion coefficient for salt and tracers (m^2/s)
80 c KPPghat - Nonlocal transport coefficient (s/m^2)
81 c KPPhbl - Boundary layer depth on "t-grid" (m)
82 c KPPfrac - Fraction of short-wave flux penetrating mixing layer
83
84 c-- KPP_CALC computes vertical viscosity and diffusivity for region
85 c (-2:sNx+3,-2:sNy+3) as required by CALC_DIFFUSIVITY and requires
86 c values of uVel, vVel, SurfaceTendencyU, SurfaceTendencyV in the
87 c region (-2:sNx+4,-2:sNy+4).
88 c Hence overlap region needs to be set OLx=4, OLy=4.
89 c When option FRUGAL_KPP is used, computation in overlap regions
90 c is replaced with exchange calls hence reducing overlap requirements
91 c to OLx=1, OLy=1.
92
93 #include "SIZE.h"
94 #include "EEPARAMS.h"
95 #include "PARAMS.h"
96 #include "DYNVARS.h"
97 #include "KPP.h"
98 #include "KPP_PARAMS.h"
99 #include "FFIELDS.h"
100 #include "GRID.h"
101
102 #ifdef ALLOW_AUTODIFF_TAMC
103 #include "tamc.h"
104 #include "tamc_keys.h"
105 INTEGER isbyte
106 PARAMETER( isbyte = 4 )
107 #else /* ALLOW_AUTODIFF_TAMC */
108 integer ikey
109 #endif /* ALLOW_AUTODIFF_TAMC */
110
111 EXTERNAL DIFFERENT_MULTIPLE
112 LOGICAL DIFFERENT_MULTIPLE
113
114 c Routine arguments
115 c bi, bj - array indices on which to apply calculations
116 c myTime - Current time in simulation
117
118 INTEGER bi, bj
119 INTEGER myThid
120 _RL myTime
121
122 #ifdef ALLOW_KPP
123
124 c Local constants
125 c minusone, p0, p5, p25, p125, p0625
126 c imin, imax, jmin, jmax - array computation indices
127
128 _RL minusone
129 parameter( minusone=-1.0)
130 _KPP_RL p0 , p5 , p25 , p125 , p0625
131 parameter( p0=0.0, p5=0.5, p25=0.25, p125=0.125, p0625=0.0625 )
132 integer imin , imax , jmin , jmax
133 #ifdef FRUGAL_KPP
134 parameter( imin=1 , imax=sNx , jmin=1 , jmax=sNy )
135 #else
136 parameter( imin=-2 , imax=sNx+3 , jmin=-2 , jmax=sNy+3 )
137 #endif
138
139 c Local arrays and variables
140 c work? (nx,ny) - horizontal working arrays
141 c ustar (nx,ny) - surface friction velocity (m/s)
142 c bo (nx,ny) - surface turbulent buoyancy forcing (m^2/s^3)
143 c bosol (nx,ny) - surface radiative buoyancy forcing (m^2/s^3)
144 c shsq (nx,ny,Nr) - local velocity shear squared
145 c at interfaces for ri_iwmix (m^2/s^2)
146 c dVsq (nx,ny,Nr) - velocity shear re surface squared
147 c at grid levels for bldepth (m^2/s^2)
148 c dbloc (nx,ny,Nr) - local delta buoyancy at interfaces
149 c for ri_iwmix and bldepth (m/s^2)
150 c Ritop (nx,ny,Nr) - numerator of bulk richardson number
151 c at grid levels for bldepth
152 c vddiff (nx,ny,Nrp2,1)- vertical viscosity on "t-grid" (m^2/s)
153 c vddiff (nx,ny,Nrp2,2)- vert. diff. on next row for temperature (m^2/s)
154 c vddiff (nx,ny,Nrp2,3)- vert. diff. on next row for salt&tracers (m^2/s)
155 c ghat (nx,ny,Nr) - nonlocal transport coefficient (s/m^2)
156 c hbl (nx,ny) - mixing layer depth (m)
157 c kmtj (nx,ny) - maximum number of wet levels in each column
158 c z0 (nx,ny) - Roughness length (m)
159 c zRef (nx,ny) - Reference depth: Hmix * epsilon (m)
160 c uRef (nx,ny) - Reference zonal velocity (m/s)
161 c vRef (nx,ny) - Reference meridional velocity (m/s)
162
163 _RL worka ( 1-OLx:sNx+OLx, 1-OLy:sNy+OLy )
164 integer work1 ( ibot:itop , jbot:jtop )
165 _KPP_RL work2 ( ibot:itop , jbot:jtop )
166 _KPP_RL ustar ( ibot:itop , jbot:jtop )
167 _KPP_RL bo ( ibot:itop , jbot:jtop )
168 _KPP_RL bosol ( ibot:itop , jbot:jtop )
169 _KPP_RL shsq ( ibot:itop , jbot:jtop , Nr )
170 _KPP_RL dVsq ( ibot:itop , jbot:jtop , Nr )
171 _KPP_RL dbloc ( ibot:itop , jbot:jtop , Nr )
172 _KPP_RL Ritop ( ibot:itop , jbot:jtop , Nr )
173 _KPP_RL vddiff( ibot:itop , jbot:jtop , 0:Nrp1, mdiff )
174 _KPP_RL ghat ( ibot:itop , jbot:jtop , Nr )
175 _KPP_RL hbl ( ibot:itop , jbot:jtop )
176 #ifdef KPP_ESTIMATE_UREF
177 _KPP_RL z0 ( ibot:itop , jbot:jtop )
178 _KPP_RL zRef ( ibot:itop , jbot:jtop )
179 _KPP_RL uRef ( ibot:itop , jbot:jtop )
180 _KPP_RL vRef ( ibot:itop , jbot:jtop )
181 #endif /* KPP_ESTIMATE_UREF */
182
183 _KPP_RL tempvar1, tempvar2
184 integer i, j, k, kp1, im1, ip1, jm1, jp1
185
186 #ifdef KPP_ESTIMATE_UREF
187 _KPP_RL dBdz1, dBdz2, ustarX, ustarY
188 #endif
189
190 c Check to see if new vertical mixing coefficient should be computed now?
191 IF ( DIFFERENT_MULTIPLE(kpp_freq,myTime,myTime-deltaTClock) .OR.
192 1 myTime .EQ. startTime ) THEN
193
194 c-----------------------------------------------------------------------
195 c prepare input arrays for subroutine "kppmix" to compute
196 c viscosity and diffusivity and ghat.
197 c All input arrays need to be in m-k-s units.
198 c
199 c note: for the computation of the bulk richardson number in the
200 c "bldepth" subroutine, gradients of velocity and buoyancy are
201 c required at every depth. in the case of very fine vertical grids
202 c (thickness of top layer < 2m), the surface reference depth must
203 c be set to zref=epsilon/2*zgrid(k), and the reference value
204 c of velocity and buoyancy must be computed as vertical average
205 c between the surface and 2*zref. in the case of coarse vertical
206 c grids zref is zgrid(1)/2., and the surface reference value is
207 c simply the surface value at zgrid(1).
208 c-----------------------------------------------------------------------
209
210 c------------------------------------------------------------------------
211 c density related quantities
212 c --------------------------
213 c
214 c work2 - density of surface layer (kg/m^3)
215 c dbloc - local buoyancy gradient at Nr interfaces
216 c g/rho{k+1,k+1} * [ drho{k,k+1}-drho{k+1,k+1} ] (m/s^2)
217 c dbsfc (stored in Ritop to conserve stack memory)
218 c - buoyancy difference with respect to the surface
219 c g * [ drho{1,k}/rho{1,k} - drho{k,k}/rho{k,k} ] (m/s^2)
220 c ttalpha (stored in vddiff(:,:,:,1) to conserve stack memory)
221 c - thermal expansion coefficient without 1/rho factor
222 c d(rho{k,k})/d(T(k)) (kg/m^3/C)
223 c ssbeta (stored in vddiff(:,:,:,2) to conserve stack memory)
224 c - salt expansion coefficient without 1/rho factor
225 c d(rho{k,k})/d(S(k)) (kg/m^3/PSU)
226 c------------------------------------------------------------------------
227
228 CALL TIMER_START('STATEKPP [KPP_CALC]', myThid)
229 CALL STATEKPP(
230 I bi, bj, myThid
231 O , work2, dbloc, Ritop
232 O , vddiff(ibot,jbot,1,1), vddiff(ibot,jbot,1,2)
233 & )
234 CALL TIMER_STOP ('STATEKPP [KPP_CALC]', myThid)
235
236 DO k = 1, Nr
237 DO j = jbot, jtop
238 DO i = ibot, itop
239 ghat(i,j,k) = dbloc(i,j,k)
240 ENDDO
241 ENDDO
242 ENDDO
243
244 #ifdef KPP_SMOOTH_DBLOC
245 c horizontally smooth dbloc with a 121 filter
246 c smooth dbloc stored in ghat to save space
247 c dbloc(k) is buoyancy gradientnote between k and k+1
248 c levels therefore k+1 mask must be used
249
250 DO k = 1, Nr-1
251 CALL KPP_SMOOTH_HORIZ (
252 I k+1, bi, bj,
253 U ghat (ibot,jbot,k) )
254 ENDDO
255
256 #endif /* KPP_SMOOTH_DBLOC */
257
258 #ifdef KPP_SMOOTH_DENS
259 c horizontally smooth density related quantities with 121 filters
260 CALL KPP_SMOOTH_HORIZ (
261 I 1, bi, bj,
262 U work2 )
263 DO k = 1, Nr
264 CALL KPP_SMOOTH_HORIZ (
265 I k+1, bi, bj,
266 U dbloc (ibot,jbot,k) )
267 CALL KPP_SMOOTH_HORIZ (
268 I k, bi, bj,
269 U Ritop (ibot,jbot,k) )
270 CALL KPP_SMOOTH_HORIZ (
271 I k, bi, bj,
272 U vddiff(ibot,jbot,k,1) )
273 CALL KPP_SMOOTH_HORIZ (
274 I k, bi, bj,
275 U vddiff(ibot,jbot,k,2) )
276 ENDDO
277 #endif /* KPP_SMOOTH_DENS */
278
279 DO k = 1, Nr
280 DO j = jbot, jtop
281 DO i = ibot, itop
282
283 c zero out dbloc over land points (so that the convective
284 c part of the interior mixing can be diagnosed)
285 dbloc(i,j,k) = dbloc(i,j,k) * pMask(i,j,k,bi,bj)
286 ghat(i,j,k) = ghat(i,j,k) * pMask(i,j,k,bi,bj)
287 Ritop(i,j,k) = Ritop(i,j,k) * pMask(i,j,k,bi,bj)
288 if(k.eq.nzmax(i,j,bi,bj)) then
289 dbloc(i,j,k) = p0
290 ghat(i,j,k) = p0
291 Ritop(i,j,k) = p0
292 endif
293
294 c numerator of bulk richardson number on grid levels
295 c note: land and ocean bottom values need to be set to zero
296 c so that the subroutine "bldepth" works correctly
297 Ritop(i,j,k) = (zgrid(1)-zgrid(k)) * Ritop(i,j,k)
298
299 END DO
300 END DO
301 END DO
302
303 c------------------------------------------------------------------------
304 c friction velocity, turbulent and radiative surface buoyancy forcing
305 c -------------------------------------------------------------------
306 c taux / rho = SurfaceTendencyU * delZ(1) (N/m^2)
307 c tauy / rho = SurfaceTendencyV * delZ(1) (N/m^2)
308 c ustar = sqrt( sqrt( taux^2 + tauy^2 ) / rho ) (m/s)
309 c bo = - g * ( alpha*SurfaceTendencyT +
310 c beta *SurfaceTendencyS ) * delZ(1) / rho (m^2/s^3)
311 c bosol = - g * alpha * Qsw * delZ(1) / rho (m^2/s^3)
312 c------------------------------------------------------------------------
313
314 c initialize arrays to zero
315 DO j = jbot, jtop
316 DO i = ibot, itop
317 ustar(i,j) = p0
318 bo (I,J) = p0
319 bosol(I,J) = p0
320 END DO
321 END DO
322
323 DO j = jmin, jmax
324 jp1 = j + 1
325 DO i = imin, imax
326 ip1 = i+1
327 tempVar1 =
328 & (SurfaceTendencyU(i,j,bi,bj) + SurfaceTendencyU(ip1,j,bi,bj)) *
329 & (SurfaceTendencyU(i,j,bi,bj) + SurfaceTendencyU(ip1,j,bi,bj)) +
330 & (SurfaceTendencyV(i,j,bi,bj) + SurfaceTendencyV(i,jp1,bi,bj)) *
331 & (SurfaceTendencyV(i,j,bi,bj) + SurfaceTendencyV(i,jp1,bi,bj))
332 if ( tempVar1 .lt. (phepsi*phepsi) ) then
333 ustar(i,j) = SQRT( phepsi * p5 * delZ(1) )
334 else
335 tempVar2 = SQRT( tempVar1 ) * p5 * delZ(1)
336 ustar(i,j) = SQRT( tempVar2 )
337 endif
338 bo(I,J) = - gravity *
339 & ( vddiff(I,J,1,1) * SurfaceTendencyT(i,j,bi,bj) +
340 & vddiff(I,J,1,2) * SurfaceTendencyS(i,j,bi,bj)
341 & ) *
342 & delZ(1) / work2(I,J)
343 bosol(I,J) = gravity * vddiff(I,J,1,1) * Qsw(i,j,bi,bj) *
344 & recip_Cp*recip_rhoNil*recip_dRf(1) *
345 & delZ(1) / work2(I,J)
346 END DO
347 END DO
348
349 c------------------------------------------------------------------------
350 c velocity shear
351 c --------------
352 c Get velocity shear squared, averaged from "u,v-grid"
353 c onto "t-grid" (in (m/s)**2):
354 c dVsq(k)=(Uref-U(k))**2+(Vref-V(k))**2 at grid levels
355 c shsq(k)=(U(k)-U(k+1))**2+(V(k)-V(k+1))**2 at interfaces
356 c------------------------------------------------------------------------
357
358 c initialize arrays to zero
359 DO k = 1, Nr
360 DO j = jbot, jtop
361 DO i = ibot, itop
362 shsq(i,j,k) = p0
363 dVsq(i,j,k) = p0
364 END DO
365 END DO
366 END DO
367
368 c dVsq computation
369
370 #ifdef KPP_ESTIMATE_UREF
371
372 c Get rid of vertical resolution dependence of dVsq term by
373 c estimating a surface velocity that is independent of first level
374 c thickness in the model. First determine mixed layer depth hMix.
375 c Second zRef = espilon * hMix. Third determine roughness length
376 c scale z0. Third estimate reference velocity.
377
378 DO j = jmin, jmax
379 jp1 = j + 1
380 DO i = imin, imax
381 ip1 = i + 1
382
383 c Determine mixed layer depth hMix as the shallowest depth at which
384 c dB/dz exceeds 5.2e-5 s^-2.
385 work1(i,j) = nzmax(i,j,bi,bj)
386 DO k = 1, Nr
387 IF ( k .LT. nzmax(i,j,bi,bj) .AND.
388 & dbloc(i,j,k) / drC(k+1) .GT. dB_dz )
389 & work1(i,j) = k
390 END DO
391
392 c Linearly interpolate to find hMix.
393 k = work1(i,j)
394 IF ( k .EQ. 0 .OR. nzmax(i,j,bi,bj) .EQ. 1 ) THEN
395 zRef(i,j) = p0
396 ELSEIF ( k .EQ. 1) THEN
397 dBdz2 = dbloc(i,j,1) / drC(2)
398 zRef(i,j) = drF(1) * dB_dz / dBdz2
399 ELSEIF ( k .LT. nzmax(i,j,bi,bj) ) THEN
400 dBdz1 = dbloc(i,j,k-1) / drC(k )
401 dBdz2 = dbloc(i,j,k ) / drC(k+1)
402 zRef(i,j) = rF(k) + drF(k) * (dB_dz - dBdz1) /
403 & MAX ( phepsi, dBdz2 - dBdz1 )
404 ELSE
405 zRef(i,j) = rF(k+1)
406 ENDIF
407
408 c Compute roughness length scale z0 subject to 0 < z0
409 tempVar1 = p5 * (
410 & (uVel(i, j, 1,bi,bj)-uVel(i, j, 2,bi,bj)) *
411 & (uVel(i, j, 1,bi,bj)-uVel(i, j, 2,bi,bj)) +
412 & (uVel(ip1,j, 1,bi,bj)-uVel(ip1,j, 2,bi,bj)) *
413 & (uVel(ip1,j, 1,bi,bj)-uVel(ip1,j, 2,bi,bj)) +
414 & (vVel(i, j, 1,bi,bj)-vVel(i, j, 2,bi,bj)) *
415 & (vVel(i, j, 1,bi,bj)-vVel(i, j, 2,bi,bj)) +
416 & (vVel(i, jp1,1,bi,bj)-vVel(i, jp1,2,bi,bj)) *
417 & (vVel(i, jp1,1,bi,bj)-vVel(i, jp1,2,bi,bj)) )
418 if ( tempVar1 .lt. (epsln*epsln) ) then
419 tempVar2 = epsln
420 else
421 tempVar2 = SQRT ( tempVar1 )
422 endif
423 z0(i,j) = rF(2) *
424 & ( rF(3) * LOG ( rF(3) / rF(2) ) /
425 & ( rF(3) - rF(2) ) -
426 & tempVar2 * vonK /
427 & MAX ( ustar(i,j), phepsi ) )
428 z0(i,j) = MAX ( z0(i,j), phepsi )
429
430 c zRef is set to 0.1 * hMix subject to z0 <= zRef <= drF(1)
431 zRef(i,j) = MAX ( epsilon * zRef(i,j), z0(i,j) )
432 zRef(i,j) = MIN ( zRef(i,j), drF(1) )
433
434 c Estimate reference velocity uRef and vRef.
435 uRef(i,j) = p5 *
436 & ( uVel(i,j,1,bi,bj) + uVel(ip1,j,1,bi,bj) )
437 vRef(i,j) = p5 *
438 & ( vVel(i,j,1,bi,bj) + vVel(i,jp1,1,bi,bj) )
439 IF ( zRef(i,j) .LT. drF(1) ) THEN
440 ustarX = ( SurfaceTendencyU(i, j,bi,bj) +
441 & SurfaceTendencyU(ip1,j,bi,bj) ) * p5
442 ustarY = ( SurfaceTendencyV(i,j, bi,bj) +
443 & SurfaceTendencyU(i,jp1,bi,bj) ) * p5
444 tempVar1 = ustarX * ustarX + ustarY * ustarY
445 if ( tempVar1 .lt. (epsln*epsln) ) then
446 tempVar2 = epsln
447 else
448 tempVar2 = SQRT ( tempVar1 )
449 endif
450 tempVar2 = ustar(i,j) *
451 & ( LOG ( zRef(i,j) / rF(2) ) +
452 & z0(i,j) / zRef(i,j) - z0(i,j) / rF(2) ) /
453 & vonK / tempVar2
454 uRef(i,j) = uRef(i,j) + ustarX * tempVar2
455 vRef(i,j) = vRef(i,j) + ustarY * tempVar2
456 ENDIF
457
458 END DO
459 END DO
460
461 DO k = 1, Nr
462 DO j = jmin, jmax
463 jm1 = j - 1
464 jp1 = j + 1
465 DO i = imin, imax
466 im1 = i - 1
467 ip1 = i + 1
468 dVsq(i,j,k) = p5 * (
469 $ (uRef(i,j) - uVel(i, j, k,bi,bj)) *
470 $ (uRef(i,j) - uVel(i, j, k,bi,bj)) +
471 $ (uRef(i,j) - uVel(ip1,j, k,bi,bj)) *
472 $ (uRef(i,j) - uVel(ip1,j, k,bi,bj)) +
473 $ (vRef(i,j) - vVel(i, j, k,bi,bj)) *
474 $ (vRef(i,j) - vVel(i, j, k,bi,bj)) +
475 $ (vRef(i,j) - vVel(i, jp1,k,bi,bj)) *
476 $ (vRef(i,j) - vVel(i, jp1,k,bi,bj)) )
477 #ifdef KPP_SMOOTH_DVSQ
478 dVsq(i,j,k) = p5 * dVsq(i,j,k) + p125 * (
479 $ (uRef(i,j) - uVel(i, jm1,k,bi,bj)) *
480 $ (uRef(i,j) - uVel(i, jm1,k,bi,bj)) +
481 $ (uRef(i,j) - uVel(ip1,jm1,k,bi,bj)) *
482 $ (uRef(i,j) - uVel(ip1,jm1,k,bi,bj)) +
483 $ (uRef(i,j) - uVel(i, jp1,k,bi,bj)) *
484 $ (uRef(i,j) - uVel(i, jp1,k,bi,bj)) +
485 $ (uRef(i,j) - uVel(ip1,jp1,k,bi,bj)) *
486 $ (uRef(i,j) - uVel(ip1,jp1,k,bi,bj)) +
487 $ (vRef(i,j) - vVel(im1,j, k,bi,bj)) *
488 $ (vRef(i,j) - vVel(im1,j, k,bi,bj)) +
489 $ (vRef(i,j) - vVel(im1,jp1,k,bi,bj)) *
490 $ (vRef(i,j) - vVel(im1,jp1,k,bi,bj)) +
491 $ (vRef(i,j) - vVel(ip1,j, k,bi,bj)) *
492 $ (vRef(i,j) - vVel(ip1,j, k,bi,bj)) +
493 $ (vRef(i,j) - vVel(ip1,jp1,k,bi,bj)) *
494 $ (vRef(i,j) - vVel(ip1,jp1,k,bi,bj)) )
495 #endif /* KPP_SMOOTH_DVSQ */
496 END DO
497 END DO
498 END DO
499
500 #else /* KPP_ESTIMATE_UREF */
501
502 DO k = 1, Nr
503 DO j = jmin, jmax
504 jm1 = j - 1
505 jp1 = j + 1
506 DO i = imin, imax
507 im1 = i - 1
508 ip1 = i + 1
509 dVsq(i,j,k) = p5 * (
510 $ (uVel(i, j, 1,bi,bj)-uVel(i, j, k,bi,bj)) *
511 $ (uVel(i, j, 1,bi,bj)-uVel(i, j, k,bi,bj)) +
512 $ (uVel(ip1,j, 1,bi,bj)-uVel(ip1,j, k,bi,bj)) *
513 $ (uVel(ip1,j, 1,bi,bj)-uVel(ip1,j, k,bi,bj)) +
514 $ (vVel(i, j, 1,bi,bj)-vVel(i, j, k,bi,bj)) *
515 $ (vVel(i, j, 1,bi,bj)-vVel(i, j, k,bi,bj)) +
516 $ (vVel(i, jp1,1,bi,bj)-vVel(i, jp1,k,bi,bj)) *
517 $ (vVel(i, jp1,1,bi,bj)-vVel(i, jp1,k,bi,bj)) )
518 #ifdef KPP_SMOOTH_DVSQ
519 dVsq(i,j,k) = p5 * dVsq(i,j,k) + p125 * (
520 $ (uVel(i, jm1,1,bi,bj)-uVel(i, jm1,k,bi,bj)) *
521 $ (uVel(i, jm1,1,bi,bj)-uVel(i, jm1,k,bi,bj)) +
522 $ (uVel(ip1,jm1,1,bi,bj)-uVel(ip1,jm1,k,bi,bj)) *
523 $ (uVel(ip1,jm1,1,bi,bj)-uVel(ip1,jm1,k,bi,bj)) +
524 $ (uVel(i, jp1,1,bi,bj)-uVel(i, jp1,k,bi,bj)) *
525 $ (uVel(i, jp1,1,bi,bj)-uVel(i, jp1,k,bi,bj)) +
526 $ (uVel(ip1,jp1,1,bi,bj)-uVel(ip1,jp1,k,bi,bj)) *
527 $ (uVel(ip1,jp1,1,bi,bj)-uVel(ip1,jp1,k,bi,bj)) +
528 $ (vVel(im1,j, 1,bi,bj)-vVel(im1,j, k,bi,bj)) *
529 $ (vVel(im1,j, 1,bi,bj)-vVel(im1,j, k,bi,bj)) +
530 $ (vVel(im1,jp1,1,bi,bj)-vVel(im1,jp1,k,bi,bj)) *
531 $ (vVel(im1,jp1,1,bi,bj)-vVel(im1,jp1,k,bi,bj)) +
532 $ (vVel(ip1,j, 1,bi,bj)-vVel(ip1,j, k,bi,bj)) *
533 $ (vVel(ip1,j, 1,bi,bj)-vVel(ip1,j, k,bi,bj)) +
534 $ (vVel(ip1,jp1,1,bi,bj)-vVel(ip1,jp1,k,bi,bj)) *
535 $ (vVel(ip1,jp1,1,bi,bj)-vVel(ip1,jp1,k,bi,bj)) )
536 #endif /* KPP_SMOOTH_DVSQ */
537 END DO
538 END DO
539 END DO
540
541 #endif /* KPP_ESTIMATE_UREF */
542
543 c shsq computation
544 DO k = 1, Nrm1
545 kp1 = k + 1
546 DO j = jmin, jmax
547 jm1 = j - 1
548 jp1 = j + 1
549 DO i = imin, imax
550 im1 = i - 1
551 ip1 = i + 1
552 shsq(i,j,k) = p5 * (
553 $ (uVel(i, j, k,bi,bj)-uVel(i, j, kp1,bi,bj)) *
554 $ (uVel(i, j, k,bi,bj)-uVel(i, j, kp1,bi,bj)) +
555 $ (uVel(ip1,j, k,bi,bj)-uVel(ip1,j, kp1,bi,bj)) *
556 $ (uVel(ip1,j, k,bi,bj)-uVel(ip1,j, kp1,bi,bj)) +
557 $ (vVel(i, j, k,bi,bj)-vVel(i, j, kp1,bi,bj)) *
558 $ (vVel(i, j, k,bi,bj)-vVel(i, j, kp1,bi,bj)) +
559 $ (vVel(i, jp1,k,bi,bj)-vVel(i, jp1,kp1,bi,bj)) *
560 $ (vVel(i, jp1,k,bi,bj)-vVel(i, jp1,kp1,bi,bj)) )
561 #ifdef KPP_SMOOTH_SHSQ
562 shsq(i,j,k) = p5 * shsq(i,j,k) + p125 * (
563 $ (uVel(i, jm1,k,bi,bj)-uVel(i, jm1,kp1,bi,bj)) *
564 $ (uVel(i, jm1,k,bi,bj)-uVel(i, jm1,kp1,bi,bj)) +
565 $ (uVel(ip1,jm1,k,bi,bj)-uVel(ip1,jm1,kp1,bi,bj)) *
566 $ (uVel(ip1,jm1,k,bi,bj)-uVel(ip1,jm1,kp1,bi,bj)) +
567 $ (uVel(i, jp1,k,bi,bj)-uVel(i, jp1,kp1,bi,bj)) *
568 $ (uVel(i, jp1,k,bi,bj)-uVel(i, jp1,kp1,bi,bj)) +
569 $ (uVel(ip1,jp1,k,bi,bj)-uVel(ip1,jp1,kp1,bi,bj)) *
570 $ (uVel(ip1,jp1,k,bi,bj)-uVel(ip1,jp1,kp1,bi,bj)) +
571 $ (vVel(im1,j, k,bi,bj)-vVel(im1,j, kp1,bi,bj)) *
572 $ (vVel(im1,j, k,bi,bj)-vVel(im1,j, kp1,bi,bj)) +
573 $ (vVel(im1,jp1,k,bi,bj)-vVel(im1,jp1,kp1,bi,bj)) *
574 $ (vVel(im1,jp1,k,bi,bj)-vVel(im1,jp1,kp1,bi,bj)) +
575 $ (vVel(ip1,j, k,bi,bj)-vVel(ip1,j, kp1,bi,bj)) *
576 $ (vVel(ip1,j, k,bi,bj)-vVel(ip1,j, kp1,bi,bj)) +
577 $ (vVel(ip1,jp1,k,bi,bj)-vVel(ip1,jp1,kp1,bi,bj)) *
578 $ (vVel(ip1,jp1,k,bi,bj)-vVel(ip1,jp1,kp1,bi,bj)) )
579 #endif
580 END DO
581 END DO
582 END DO
583
584 c-----------------------------------------------------------------------
585 c solve for viscosity, diffusivity, ghat, and hbl on "t-grid"
586 c-----------------------------------------------------------------------
587
588 DO j = jbot, jtop
589 DO i = ibot, itop
590 work1(i,j) = nzmax(i,j,bi,bj)
591 work2(i,j) = Fcori(i,j,bi,bj)
592 END DO
593 END DO
594 CALL TIMER_START('KPPMIX [KPP_CALC]', myThid)
595 CALL KPPMIX (
596 I mytime, mythid
597 I , work1, shsq, dVsq, ustar
598 I , bo, bosol, dbloc, Ritop, work2
599 I , ikey
600 O , vddiff
601 U , ghat
602 O , hbl )
603
604 CALL TIMER_STOP ('KPPMIX [KPP_CALC]', myThid)
605
606 #ifdef ALLOW_AUTODIFF_TAMC
607 cph( storing not necessary
608 cphCADJ STORE vddiff, ghat = comlev1_kpp, key = ikey
609 cph)
610 #endif /* ALLOW_AUTODIFF_TAMC */
611
612 c-----------------------------------------------------------------------
613 c zero out land values and transfer to global variables
614 c-----------------------------------------------------------------------
615
616 DO j = jmin, jmax
617 DO i = imin, imax
618 DO k = 1, Nr
619 KPPviscAz(i,j,k,bi,bj) = vddiff(i,j,k-1,1) * pMask(i,j,k,bi,bj)
620 KPPdiffKzS(i,j,k,bi,bj)= vddiff(i,j,k-1,2) * pMask(i,j,k,bi,bj)
621 KPPdiffKzT(i,j,k,bi,bj)= vddiff(i,j,k-1,3) * pMask(i,j,k,bi,bj)
622 KPPghat(i,j,k,bi,bj) = ghat(i,j,k) * pMask(i,j,k,bi,bj)
623 END DO
624 KPPhbl(i,j,bi,bj) = hbl(i,j) * pMask(i,j,1,bi,bj)
625 END DO
626 END DO
627 #ifdef FRUGAL_KPP
628 _EXCH_XYZ_R8(KPPviscAz , myThid )
629 _EXCH_XYZ_R8(KPPdiffKzS , myThid )
630 _EXCH_XYZ_R8(KPPdiffKzT , myThid )
631 _EXCH_XYZ_R8(KPPghat , myThid )
632 _EXCH_XY_R8 (KPPhbl , myThid )
633 #endif
634
635 #ifdef KPP_SMOOTH_VISC
636 c horizontal smoothing of vertical viscosity
637 DO k = 1, Nr
638 CALL SMOOTH_HORIZ (
639 I k, bi, bj,
640 U KPPviscAz(1-OLx,1-OLy,k,bi,bj) )
641 END DO
642 _EXCH_XYZ_R8(KPPviscAz , myThid )
643 #endif /* KPP_SMOOTH_VISC */
644
645 #ifdef KPP_SMOOTH_DIFF
646 c horizontal smoothing of vertical diffusivity
647 DO k = 1, Nr
648 CALL SMOOTH_HORIZ (
649 I k, bi, bj,
650 U KPPdiffKzS(1-OLx,1-OLy,k,bi,bj) )
651 CALL SMOOTH_HORIZ (
652 I k, bi, bj,
653 U KPPdiffKzT(1-OLx,1-OLy,k,bi,bj) )
654 END DO
655 _EXCH_XYZ_R8(KPPdiffKzS , myThid )
656 _EXCH_XYZ_R8(KPPdiffKzT , myThid )
657 #endif /* KPP_SMOOTH_DIFF */
658
659 C Compute fraction of solar short-wave flux penetrating to
660 C the bottom of the mixing layer.
661 DO j=1-OLy,sNy+OLy
662 DO i=1-OLx,sNx+OLx
663 worka(i,j) = KPPhbl(i,j,bi,bj)
664 ENDDO
665 ENDDO
666 CALL SWFRAC(
667 I (sNx+2*OLx)*(sNy+2*OLy), minusone,
668 I mytime, mythid,
669 U worka )
670 DO j=1-OLy,sNy+OLy
671 DO i=1-OLx,sNx+OLx
672 KPPfrac(i,j,bi,bj) = worka(i,j)
673 ENDDO
674 ENDDO
675
676 ENDIF
677
678 #endif ALLOW_KPP
679
680 RETURN
681 END

  ViewVC Help
Powered by ViewVC 1.1.22