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

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

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


Revision 1.22 - (hide annotations) (download)
Thu Nov 6 22:12:09 2003 UTC (20 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint52, checkpoint52a_pre, ecco_c52_e35, checkpoint51u_post
Changes since 1.21: +8 -5 lines
o merging from ecco-branch

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

  ViewVC Help
Powered by ViewVC 1.1.22