/[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.4 - (hide annotations) (download)
Mon Nov 13 16:37:02 2000 UTC (24 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint32
Changes since 1.3: +132 -334 lines
Modified and fixed version. Tested with verification/natl_box.

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

  ViewVC Help
Powered by ViewVC 1.1.22