/[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.29 - (show annotations) (download)
Sun Jul 18 15:29:52 2004 UTC (19 years, 10 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint54d_post, checkpoint54e_post, checkpoint55, checkpoint54f_post, checkpoint55c_post, checkpoint55g_post, checkpoint55d_post, checkpoint55d_pre, checkpoint55b_post, checkpoint55f_post, checkpoint55a_post, checkpoint55e_post
Changes since 1.28: +2 -2 lines
bug fixed in kpp_calc.F : ustarY = ( SurfaceTendencyV )_j /drF(1)

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

  ViewVC Help
Powered by ViewVC 1.1.22