/[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.23 - (show annotations) (download)
Thu Nov 13 06:35:14 2003 UTC (20 years, 6 months ago) by dimitri
Branch: MAIN
CVS Tags: checkpoint52l_pre, hrcube4, hrcube5, checkpoint52d_pre, checkpoint52j_pre, checkpoint52l_post, checkpoint52k_post, checkpoint52f_post, checkpoint52i_pre, hrcube_1, hrcube_2, hrcube_3, checkpoint52e_pre, checkpoint52e_post, checkpoint52b_pre, checkpoint52m_post, checkpoint52b_post, checkpoint52c_post, checkpoint52f_pre, checkpoint52d_post, checkpoint52i_post, checkpoint52h_pre, checkpoint52j_post, branch-netcdf, checkpoint52a_post
Branch point for: netcdf-sm0
Changes since 1.22: +4 -22 lines
o modifications to make FREEZE flux visible to pkg/kpp
  - moved surfaceTendencyTice from pkg/seaice to main code
  - FREEZE moved to FORWARD_STEP
  - subroutine FREEZE now limits only surface temperature
    this means new output.txt for global_ocean.90x40x15,
    global_ocean.cs32x15, and global_with_exf, but note
    that results for these three experiments remain
    bit-identical to before if allowFreezing=.FALSE.
o added surface flux output variables to TIMEAVE_STATVARS
o time-averaged output for pkg/ptracers

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

  ViewVC Help
Powered by ViewVC 1.1.22