/[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.34 - (show annotations) (download)
Thu Apr 19 04:51:59 2007 UTC (18 years, 3 months ago) by dimitri
Branch: MAIN
CVS Tags: checkpoint59
Changes since 1.33: +35 -49 lines
FRUGAL_KPP option removed due to popular undemand

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

  ViewVC Help
Powered by ViewVC 1.1.22