/[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.3 - (hide annotations) (download)
Wed Sep 13 17:07:11 2000 UTC (23 years, 8 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint31
Changes since 1.2: +0 -0 lines
KPP package without flag kppPackageIsOn.

1 heimbach 1.2 C $Header: /escher1/cvs/master/mitgcmuv/pkg/kpp/kpp_calc.F,v 1.5 2000/09/11 22:32:53 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     c Local arrays and variables
125     c work? (nx,ny) - horizontal working arrays
126     c ustar (nx,ny) - surface friction velocity (m/s)
127     c bo (nx,ny) - surface turbulent buoyancy forcing (m^2/s^3)
128     c bosol (nx,ny) - surface radiative buoyancy forcing (m^2/s^3)
129     c shsq (nx,ny,Nr) - local velocity shear squared
130     c at interfaces for ri_iwmix (m^2/s^2)
131     c dVsq (nx,ny,Nr) - velocity shear re surface squared
132     c at grid levels for bldepth (m^2/s^2)
133     c dbloc (nx,ny,Nr) - local delta buoyancy at interfaces
134     c for ri_iwmix and bldepth (m/s^2)
135     c Ritop (nx,ny,Nr) - numerator of bulk richardson number
136     c at grid levels for bldepth
137     c vddiff (nx,ny,Nrp2,1)- vertical viscosity on "t-grid" (m^2/s)
138     c vddiff (nx,ny,Nrp2,2)- vert. diff. on next row for temperature (m^2/s)
139     c vddiff (nx,ny,Nrp2,3)- vert. diff. on next row for salt&tracers (m^2/s)
140     c ghat (nx,ny,Nr) - nonlocal transport coefficient (s/m^2)
141     c hbl (nx,ny) - mixing layer depth (m)
142     c kmtj (nx,ny) - maximum number of wet levels in each column
143     c z0 (nx,ny) - Roughness length (m)
144     c zRef (nx,ny) - Reference depth: Hmix * epsilon (m)
145     c uRef (nx,ny) - Reference zonal velocity (m/s)
146     c vRef (nx,ny) - Reference meridional velocity (m/s)
147    
148     _RS worka (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
149     _RS workb (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
150     #ifdef FRUGAL_KPP
151     integer work1(sNx,sNy)
152     _RS work2 (sNx,sNy)
153     _RS ustar (sNx,sNy)
154     _RS bo (sNx,sNy)
155     _RS bosol (sNx,sNy)
156     _RS shsq (sNx,sNy,Nr)
157     _RS dVsq (sNx,sNy,Nr)
158     _RS dbloc (sNx,sNy,Nr)
159     _RS Ritop (sNx,sNy,Nr)
160     _RS vddiff (sNx,sNy,0:Nrp1,mdiff)
161     _RS ghat (sNx,sNy,Nr)
162     _RS hbl (sNx,sNy)
163     #ifdef KPP_ESTIMATE_UREF
164     _RS z0 (sNx,sNy)
165     _RS zRef (sNx,sNy)
166     _RS uRef (sNx,sNy)
167     _RS vRef (sNx,sNy)
168     #endif /* KPP_ESTIMATE_UREF */
169     #else /* FRUGAL_KPP */
170     integer work1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
171     _RS work2 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
172     _RS ustar (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
173     _RS bo (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
174     _RS bosol (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
175     _RS shsq (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
176     _RS dVsq (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
177     _RS dbloc (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
178     _RS Ritop (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
179     _RS vddiff (1-OLx:sNx+OLx,1-OLy:sNy+OLy,0:Nrp1,mdiff)
180     _RS ghat (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
181     _RS hbl (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
182     #ifdef KPP_ESTIMATE_UREF
183     _RS z0 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
184     _RS zRef (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
185     _RS uRef (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
186     _RS vRef (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
187     #endif /* KPP_ESTIMATE_UREF */
188     #endif /* FRUGAL_KPP */
189    
190     c imin,imax,jmin,jmax - array indices
191     integer imin , imax , jmin , jmax
192     parameter( imin=-2, imax=sNx+3, jmin=-2, jmax=sNy+3 )
193    
194     c mixing process switches
195     logical lri
196     parameter( lri = .true. )
197    
198 heimbach 1.2 _RS m1
199     parameter( m1=-1.0)
200 adcroft 1.1 _RS p0 , p5 , p25 , p125 , p0625
201     parameter( p0=0.0, p5=0.5, p25=0.25, p125=0.125, p0625=0.0625 )
202    
203 heimbach 1.2 _RL tempVar
204 adcroft 1.1 integer i, j, k, kp1, im1, ip1, jm1, jp1
205    
206     #ifdef KPP_ESTIMATE_UREF
207     _RS dBdz1, dBdz2, ustarX, ustarY
208     #endif
209    
210     c Check to see if new vertical mixing coefficient should be computed now?
211     IF ( DIFFERENT_MULTIPLE(kpp_freq,myTime,myTime-deltaTClock) .OR.
212     1 myTime .EQ. startTime ) THEN
213    
214     c-----------------------------------------------------------------------
215     c prepare input arrays for subroutine "kppmix" to compute
216     c viscosity and diffusivity and ghat.
217     c All input arrays need to be in m-k-s units.
218     c
219     c note: for the computation of the bulk richardson number in the
220     c "bldepth" subroutine, gradients of velocity and buoyancy are
221     c required at every depth. in the case of very fine vertical grids
222     c (thickness of top layer < 2m), the surface reference depth must
223     c be set to zref=epsilon/2*zgrid(k), and the reference value
224     c of velocity and buoyancy must be computed as vertical average
225     c between the surface and 2*zref. in the case of coarse vertical
226     c grids zref is zgrid(1)/2., and the surface reference value is
227     c simply the surface value at zgrid(1).
228     c-----------------------------------------------------------------------
229    
230     c------------------------------------------------------------------------
231     c density related quantities
232     c --------------------------
233     c
234     c work2 - density of surface layer (kg/m^3)
235     c dbloc - local buoyancy gradient at Nr interfaces
236     c g/rho{k+1,k+1} * [ drho{k,k+1}-drho{k+1,k+1} ] (m/s^2)
237     c dbsfc (stored in Ritop to conserve stack memory)
238     c - buoyancy difference with respect to the surface
239     c g * [ drho{1,k}/rho{1,k} - drho{k,k}/rho{k,k} ] (m/s^2)
240     c ttalpha (stored in vddiff(:,:,:,1) to conserve stack memory)
241     c - thermal expansion coefficient without 1/rho factor
242     c d(rho{k,k})/d(T(k)) (kg/m^3/C)
243     c ssbeta (stored in vddiff(:,:,:,2) to conserve stack memory)
244     c - salt expansion coefficient without 1/rho factor
245     c d(rho{k,k})/d(S(k)) (kg/m^3/PSU)
246     c------------------------------------------------------------------------
247    
248     CALL TIMER_START('STATEKPP [KPP_CALC]', myThid)
249     CALL STATEKPP(
250     I bi, bj, myThid
251     O , work2, dbloc, Ritop
252     #ifdef FRUGAL_KPP
253     O , vddiff(1 ,1 ,1,1), vddiff(1 ,1 ,1,2)
254     #else
255     O , vddiff(1-OLx,1-OLy,1,1), vddiff(1-OLx,1-OLy,1,2)
256     #endif
257     & )
258     CALL TIMER_STOP ('STATEKPP [KPP_CALC]', myThid)
259    
260     #ifdef KPP_SMOOTH_DBLOC
261     c horizontally smooth dbloc with a 121 filter
262     c (stored in ghat to save space)
263    
264     DO k = 1, Nr
265     CALL SMOOTH_HORIZ_RS (
266     I k, bi, bj,
267     I dbloc(1-OLx,1-OLy,k),
268     O ghat (1-OLx,1-OLy,k) )
269     ENDDO
270    
271     #else /* KPP_SMOOTH_DBLOC */
272    
273     DO k = 1, Nr
274     #ifdef FRUGAL_KPP
275     DO j = 1, sNy
276     DO i = 1, sNx
277     #else
278     DO j = 1-OLy, sNy+OLy
279     DO i = imin, imax
280     #endif
281     ghat(i,j,k) = dbloc(i,j,k)
282     ENDDO
283     ENDDO
284     ENDDO
285    
286     #endif /* KPP_SMOOTH_DBLOC */
287    
288     #ifdef KPP_SMOOTH_DENS
289     c horizontally smooth density related quantities with 121 filters
290     CALL SMOOTH_HORIZ_RS (
291     I k, bi, bj,
292     I work2,
293     O work2 )
294     DO k = 1, Nr
295     CALL SMOOTH_HORIZ_RS (
296     I k, bi, bj,
297     I dbloc (1-OLx,1-OLy,k) ,
298     O dbloc (1-OLx,1-OLy,k) )
299     CALL SMOOTH_HORIZ_RS (
300     I k, bi, bj,
301     I Ritop (1-OLx,1-OLy,k) ,
302     O Ritop (1-OLx,1-OLy,k) )
303     CALL SMOOTH_HORIZ_RS (
304     I k, bi, bj,
305     I vddiff(1-OLx,1-OLy,k,1),
306     O vddiff(1-OLx,1-OLy,k,1) )
307     CALL SMOOTH_HORIZ_RS (
308     I k, bi, bj,
309     I vddiff(1-OLx,1-OLy,k,2),
310     O vddiff(1-OLx,1-OLy,k,2) )
311     ENDDO
312     #endif /* KPP_SMOOTH_DENS */
313    
314     DO k = 1, Nr
315     #ifdef FRUGAL_KPP
316     DO j = 1, sNy
317     DO i = 1, sNx
318     #else
319     DO j = 1-OLy, sNy+OLy
320     DO i = 1-OLx, sNx+OLx
321     #endif
322    
323     c zero out dbloc over land points (so that the convective
324     c part of the interior mixing can be diagnosed)
325     dbloc(i,j,k) = dbloc(i,j,k) * pMask(i,j,k,bi,bj)
326     ghat(i,j,k) = ghat(i,j,k) * pMask(i,j,k,bi,bj)
327     Ritop(i,j,k) = Ritop(i,j,k) * pMask(i,j,k,bi,bj)
328     if(k.eq.nzmax(i,j,bi,bj)) then
329     dbloc(i,j,k) = p0
330     ghat(i,j,k) = p0
331     Ritop(i,j,k) = p0
332     endif
333    
334     c numerator of bulk richardson number on grid levels
335     c note: land and ocean bottom values need to be set to zero
336     c so that the subroutine "bldepth" works correctly
337     Ritop(i,j,k) = (zgrid(1)-zgrid(k)) * Ritop(i,j,k)
338    
339     END DO
340     END DO
341     END DO
342    
343     c------------------------------------------------------------------------
344     c friction velocity, turbulent and radiative surface buoyancy forcing
345     c -------------------------------------------------------------------
346 heimbach 1.2 c taux / rho = SurfaceTendencyU * delZ(1) (N/m^2)
347     c tauy / rho = SurfaceTendencyV * delZ(1) (N/m^2)
348 adcroft 1.1 c ustar = sqrt( sqrt( taux^2 + tauy^2 ) / rho ) (m/s)
349 heimbach 1.2 c bo = - g * ( alpha*SurfaceTendencyT +
350     c beta *SurfaceTendencyS ) * delZ(1) / rho (m^2/s^3)
351 adcroft 1.1 c bosol = - g * alpha * Qsw * delZ(1) / rho (m^2/s^3)
352     c------------------------------------------------------------------------
353    
354     #ifdef FRUGAL_KPP
355     DO j = 1, sNy
356     jp1 = j + 1
357     DO i = 1, sNx
358     #else
359     DO j = jmin, jmax
360 heimbach 1.2 jp1 = j + 1
361     DO i = imin, imax
362 adcroft 1.1 #endif
363 heimbach 1.2 ip1 = i+1
364     tempVar =
365     & (SurfaceTendencyU(i,j,bi,bj) + SurfaceTendencyU(ip1,j,bi,bj)) *
366     & (SurfaceTendencyU(i,j,bi,bj) + SurfaceTendencyU(ip1,j,bi,bj)) +
367     & (SurfaceTendencyV(i,j,bi,bj) + SurfaceTendencyV(i,jp1,bi,bj)) *
368     & (SurfaceTendencyV(i,j,bi,bj) + SurfaceTendencyV(i,jp1,bi,bj))
369     if ( tempVar .lt. (epsln*epsln) ) then
370     ustar(i,j) = SQRT( epsln * p5 * delZ(1) )
371     else
372     ustar(i,j) = SQRT( SQRT( tempVar ) * p5 * delZ(1) )
373     endif
374     bo(I,J) = - gravity *
375     & ( vddiff(I,J,1,1) * SurfaceTendencyT(i,j,bi,bj) +
376     & vddiff(I,J,1,2) * SurfaceTendencyS(i,j,bi,bj)
377     & ) *
378     & delZ(1) / work2(I,J)
379     bosol(I,J) = - gravity * vddiff(I,J,1,1) * Qsw(i,j,bi,bj) *
380     & delZ(1) / work2(I,J)
381     END DO
382 adcroft 1.1 END DO
383    
384     #ifndef FRUGAL_KPP
385     c set array edges to zero
386     DO j = jmin, jmax
387     DO i = 1-OLx, imin-1
388     ustar(i,j) = p0
389     bo (I,J) = p0
390     bosol(I,J) = p0
391     END DO
392     DO i = imax+1, sNx+OLx
393     ustar(i,j) = p0
394     bo (I,J) = p0
395     bosol(I,J) = p0
396     END DO
397     END DO
398     DO i = 1-OLx, sNx+OLx
399     DO j = 1-OLy, jmin-1
400     ustar(i,j) = p0
401     bo (I,J) = p0
402     bosol(I,J) = p0
403     END DO
404     DO j = jmax+1, sNy+OLy
405     ustar(i,j) = p0
406     bo (I,J) = p0
407     bosol(I,J) = p0
408     END DO
409     END DO
410     #endif
411    
412     c------------------------------------------------------------------------
413     c velocity shear
414     c --------------
415     c Get velocity shear squared, averaged from "u,v-grid"
416     c onto "t-grid" (in (m/s)**2):
417     c dVsq(k)=(Uref-U(k))**2+(Vref-V(k))**2 at grid levels
418     c shsq(k)=(U(k)-U(k+1))**2+(V(k)-V(k+1))**2 at interfaces
419     c------------------------------------------------------------------------
420    
421     c dVsq computation
422    
423     #ifdef KPP_ESTIMATE_UREF
424    
425     c Get rid of vertical resolution dependence of dVsq term by
426     c estimating a surface velocity that is independent of first level
427     c thickness in the model. First determine mixed layer depth hMix.
428     c Second zRef = espilon * hMix. Third determine roughness length
429     c scale z0. Third estimate reference velocity.
430    
431     #ifdef FRUGAL_KPP
432     DO j = 1, sNy
433     jp1 = j + 1
434     DO i = 1, sNx
435     #else
436     DO j = jmin, jmax
437     jp1 = j + 1
438     DO i = imin, imax
439     #endif /* FRUGAL_KPP */
440     ip1 = i + 1
441    
442     c Determine mixed layer depth hMix as the shallowest depth at which
443     c dB/dz exceeds 5.2e-5 s^-2.
444     work1(i,j) = nzmax(i,j,bi,bj)
445     DO k = 1, Nr
446     IF ( k .LT. nzmax(i,j,bi,bj) .AND.
447     & dbloc(i,j,k) / drC(k+1) .GT. dB_dz )
448     & work1(i,j) = k
449     END DO
450    
451     c Linearly interpolate to find hMix.
452     k = work1(i,j)
453     IF ( k .EQ. 0 .OR. nzmax(i,j,bi,bj) .EQ. 1 ) THEN
454     zRef(i,j) = p0
455     ELSEIF ( k .EQ. 1) THEN
456     dBdz2 = dbloc(i,j,1) / drC(2)
457     zRef(i,j) = drF(1) * dB_dz / dBdz2
458     ELSEIF ( k .LT. nzmax(i,j,bi,bj) ) THEN
459     dBdz1 = dbloc(i,j,k-1) / drC(k )
460     dBdz2 = dbloc(i,j,k ) / drC(k+1)
461     zRef(i,j) = rF(k) + drF(k) * (dB_dz - dBdz1) /
462     & MAX ( phepsi, dBdz2 - dBdz1 )
463     ELSE
464     zRef(i,j) = rF(k+1)
465     ENDIF
466    
467     c Compute roughness length scale z0 subject to 0 < z0
468 heimbach 1.2 tempVar = p5 * (
469 adcroft 1.1 & (uVel(i, j, 1,bi,bj)-uVel(i, j, 2,bi,bj)) *
470     & (uVel(i, j, 1,bi,bj)-uVel(i, j, 2,bi,bj)) +
471     & (uVel(ip1,j, 1,bi,bj)-uVel(ip1,j, 2,bi,bj)) *
472     & (uVel(ip1,j, 1,bi,bj)-uVel(ip1,j, 2,bi,bj)) +
473     & (vVel(i, j, 1,bi,bj)-vVel(i, j, 2,bi,bj)) *
474     & (vVel(i, j, 1,bi,bj)-vVel(i, j, 2,bi,bj)) +
475     & (vVel(i, jp1,1,bi,bj)-vVel(i, jp1,2,bi,bj)) *
476 heimbach 1.2 & (vVel(i, jp1,1,bi,bj)-vVel(i, jp1,2,bi,bj)) )
477     if ( tempVar .lt. (epsln*epsln) ) then
478     tempVar = epsln
479     else
480     tempVar = SQRT ( tempVar )
481     endif
482 adcroft 1.1 z0(i,j) = rF(2) *
483     & ( rF(3) * LOG ( rF(3) / rF(2) ) /
484     & ( rF(3) - rF(2) ) -
485     & tempVar * vonK /
486     & MAX ( ustar(i,j), phepsi ) )
487     z0(i,j) = MAX ( z0(i,j), phepsi )
488    
489     c zRef is set to 0.1 * hMix subject to z0 <= zRef <= drF(1)
490     zRef(i,j) = MAX ( epsilon * zRef(i,j), z0(i,j) )
491     zRef(i,j) = MIN ( zRef(i,j), drF(1) )
492    
493     c Estimate reference velocity uRef and vRef.
494     uRef(i,j) = p5 *
495     & ( uVel(i,j,1,bi,bj) + uVel(ip1,j,1,bi,bj) )
496     vRef(i,j) = p5 *
497     & ( vVel(i,j,1,bi,bj) + vVel(i,jp1,1,bi,bj) )
498     IF ( zRef(i,j) .LT. drF(1) ) THEN
499 heimbach 1.2 ustarX = ( SurfaceTendencyU(i, j,bi,bj) +
500     & SurfaceTendencyU(ip1,j,bi,bj) ) * p5
501     ustarY = ( SurfaceTendencyV(i,j, bi,bj) +
502     & SurfaceTendencyU(i,jp1,bi,bj) ) * p5
503     tempVar = ustarX * ustarX + ustarY * ustarY
504     if ( tempVar .lt. (epsln*epsln) ) then
505     tempVar = epsln
506     else
507     tempVar = SQRT ( tempVar )
508     endif
509 adcroft 1.1 tempVar = ustar(i,j) *
510     & ( LOG ( zRef(i,j) / rF(2) ) +
511     & z0(i,j) / zRef(i,j) - z0(i,j) / rF(2) ) /
512     & vonK / tempVar
513     uRef(i,j) = uRef(i,j) + ustarX * tempVar
514     vRef(i,j) = vRef(i,j) + ustarY * tempVar
515     ENDIF
516    
517     END DO
518     END DO
519    
520     IF (KPPmixingMaps) THEN
521     #ifdef FRUGAL_KPP
522     CALL PRINT_MAPRS(
523     I zRef, 'zRef', PRINT_MAP_XY,
524     I 1, sNx, 1, sNy, 1, 1, 1, 1,
525     I 1, sNx, 1, sNy, 1, -1, 1, 1, 1, 1, 1, 1, 1, 1, 1 )
526     CALL PRINT_MAPRS(
527     I z0, 'z0', PRINT_MAP_XY,
528     I 1, sNx, 1, sNy, 1, 1, 1, 1,
529     I 1, sNx, 1, sNy, 1, -1, 1, 1, 1, 1, 1, 1, 1, 1, 1 )
530     CALL PRINT_MAPRS(
531     I uRef, 'uRef', PRINT_MAP_XY,
532     I 1, sNx, 1, sNy, 1, 1, 1, 1,
533     I 1, sNx, 1, sNy, 1, -1, 1, 1, 1, 1, 1, 1, 1, 1, 1 )
534     CALL PRINT_MAPRS(
535     I vRef, 'vRef', PRINT_MAP_XY,
536     I 1, sNx, 1, sNy, 1, 1, 1, 1,
537     I 1, sNx, 1, sNy, 1, -1, 1, 1, 1, 1, 1, 1, 1, 1, 1 )
538     #else
539     CALL PRINT_MAPRS(
540     I zRef, 'zRef', PRINT_MAP_XY,
541     I 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1, 1, 1, 1,
542     I 1, sNx, 1, sNy, 1, -1, 1, 1, 1, 1, 1, 1, 1, 1, 1 )
543     CALL PRINT_MAPRS(
544     I z0, 'z0', PRINT_MAP_XY,
545     I 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1, 1, 1, 1,
546     I 1, sNx, 1, sNy, 1, -1, 1, 1, 1, 1, 1, 1, 1, 1, 1 )
547     CALL PRINT_MAPRS(
548     I uRef, 'uRef', PRINT_MAP_XY,
549     I 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1, 1, 1, 1,
550     I 1, sNx, 1, sNy, 1, -1, 1, 1, 1, 1, 1, 1, 1, 1, 1 )
551     CALL PRINT_MAPRS(
552     I vRef, 'vRef', PRINT_MAP_XY,
553     I 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1, 1, 1, 1,
554     I 1, sNx, 1, sNy, 1, -1, 1, 1, 1, 1, 1, 1, 1, 1, 1 )
555     #endif
556     ENDIF
557    
558     DO k = 1, Nr
559     #ifdef FRUGAL_KPP
560     DO j = 1, sNy
561     jm1 = j - 1
562     jp1 = j + 1
563     DO i = 1, sNx
564     #else
565     DO j = jmin, jmax
566     jm1 = j - 1
567     jp1 = j + 1
568     DO i = imin, imax
569     #endif /* FRUGAL_KPP */
570     im1 = i - 1
571     ip1 = i + 1
572     dVsq(i,j,k) = p5 * (
573     $ (uRef(i,j) - uVel(i, j, k,bi,bj)) *
574     $ (uRef(i,j) - uVel(i, j, k,bi,bj)) +
575     $ (uRef(i,j) - uVel(ip1,j, k,bi,bj)) *
576     $ (uRef(i,j) - uVel(ip1,j, k,bi,bj)) +
577     $ (vRef(i,j) - vVel(i, j, k,bi,bj)) *
578     $ (vRef(i,j) - vVel(i, j, k,bi,bj)) +
579     $ (vRef(i,j) - vVel(i, jp1,k,bi,bj)) *
580     $ (vRef(i,j) - vVel(i, jp1,k,bi,bj)) )
581     #ifdef KPP_SMOOTH_DVSQ
582     dVsq(i,j,k) = p5 * dVsq(i,j,k) + p125 * (
583     $ (uRef(i,j) - uVel(i, jm1,k,bi,bj)) *
584     $ (uRef(i,j) - uVel(i, jm1,k,bi,bj)) +
585     $ (uRef(i,j) - uVel(ip1,jm1,k,bi,bj)) *
586     $ (uRef(i,j) - uVel(ip1,jm1,k,bi,bj)) +
587     $ (uRef(i,j) - uVel(i, jp1,k,bi,bj)) *
588     $ (uRef(i,j) - uVel(i, jp1,k,bi,bj)) +
589     $ (uRef(i,j) - uVel(ip1,jp1,k,bi,bj)) *
590     $ (uRef(i,j) - uVel(ip1,jp1,k,bi,bj)) +
591     $ (vRef(i,j) - vVel(im1,j, k,bi,bj)) *
592     $ (vRef(i,j) - vVel(im1,j, k,bi,bj)) +
593     $ (vRef(i,j) - vVel(im1,jp1,k,bi,bj)) *
594     $ (vRef(i,j) - vVel(im1,jp1,k,bi,bj)) +
595     $ (vRef(i,j) - vVel(ip1,j, k,bi,bj)) *
596     $ (vRef(i,j) - vVel(ip1,j, k,bi,bj)) +
597     $ (vRef(i,j) - vVel(ip1,jp1,k,bi,bj)) *
598     $ (vRef(i,j) - vVel(ip1,jp1,k,bi,bj)) )
599     #endif /* KPP_SMOOTH_DVSQ */
600     END DO
601     END DO
602     END DO
603    
604     #else /* KPP_ESTIMATE_UREF */
605    
606     DO k = 1, Nr
607     #ifdef FRUGAL_KPP
608     DO j = 1, sNy
609     jm1 = j - 1
610     jp1 = j + 1
611     DO i = 1, sNx
612     #else
613     DO j = jmin, jmax
614     jm1 = j - 1
615     jp1 = j + 1
616     DO i = imin, imax
617     #endif /* FRUGAL_KPP */
618     im1 = i - 1
619     ip1 = i + 1
620     dVsq(i,j,k) = p5 * (
621     $ (uVel(i, j, 1,bi,bj)-uVel(i, j, k,bi,bj)) *
622     $ (uVel(i, j, 1,bi,bj)-uVel(i, j, k,bi,bj)) +
623     $ (uVel(ip1,j, 1,bi,bj)-uVel(ip1,j, k,bi,bj)) *
624     $ (uVel(ip1,j, 1,bi,bj)-uVel(ip1,j, k,bi,bj)) +
625     $ (vVel(i, j, 1,bi,bj)-vVel(i, j, k,bi,bj)) *
626     $ (vVel(i, j, 1,bi,bj)-vVel(i, j, k,bi,bj)) +
627     $ (vVel(i, jp1,1,bi,bj)-vVel(i, jp1,k,bi,bj)) *
628     $ (vVel(i, jp1,1,bi,bj)-vVel(i, jp1,k,bi,bj)) )
629     #ifdef KPP_SMOOTH_DVSQ
630     dVsq(i,j,k) = p5 * dVsq(i,j,k) + p125 * (
631     $ (uVel(i, jm1,1,bi,bj)-uVel(i, jm1,k,bi,bj)) *
632     $ (uVel(i, jm1,1,bi,bj)-uVel(i, jm1,k,bi,bj)) +
633     $ (uVel(ip1,jm1,1,bi,bj)-uVel(ip1,jm1,k,bi,bj)) *
634     $ (uVel(ip1,jm1,1,bi,bj)-uVel(ip1,jm1,k,bi,bj)) +
635     $ (uVel(i, jp1,1,bi,bj)-uVel(i, jp1,k,bi,bj)) *
636     $ (uVel(i, jp1,1,bi,bj)-uVel(i, jp1,k,bi,bj)) +
637     $ (uVel(ip1,jp1,1,bi,bj)-uVel(ip1,jp1,k,bi,bj)) *
638     $ (uVel(ip1,jp1,1,bi,bj)-uVel(ip1,jp1,k,bi,bj)) +
639     $ (vVel(im1,j, 1,bi,bj)-vVel(im1,j, k,bi,bj)) *
640     $ (vVel(im1,j, 1,bi,bj)-vVel(im1,j, k,bi,bj)) +
641     $ (vVel(im1,jp1,1,bi,bj)-vVel(im1,jp1,k,bi,bj)) *
642     $ (vVel(im1,jp1,1,bi,bj)-vVel(im1,jp1,k,bi,bj)) +
643     $ (vVel(ip1,j, 1,bi,bj)-vVel(ip1,j, k,bi,bj)) *
644     $ (vVel(ip1,j, 1,bi,bj)-vVel(ip1,j, k,bi,bj)) +
645     $ (vVel(ip1,jp1,1,bi,bj)-vVel(ip1,jp1,k,bi,bj)) *
646     $ (vVel(ip1,jp1,1,bi,bj)-vVel(ip1,jp1,k,bi,bj)) )
647     #endif /* KPP_SMOOTH_DVSQ */
648     END DO
649     END DO
650     END DO
651    
652     #endif /* KPP_ESTIMATE_UREF */
653    
654     c shsq computation
655     DO k = 1, Nrm1
656     kp1 = k + 1
657     #ifdef FRUGAL_KPP
658     DO j = 1, sNy
659     jm1 = j - 1
660     jp1 = j + 1
661     DO i = 1, sNx
662     #else
663     DO j = jmin, jmax
664     jm1 = j - 1
665     jp1 = j + 1
666     DO i = imin, imax
667     #endif /* FRUGAL_KPP */
668     im1 = i - 1
669     ip1 = i + 1
670     shsq(i,j,k) = p5 * (
671     $ (uVel(i, j, k,bi,bj)-uVel(i, j, kp1,bi,bj)) *
672     $ (uVel(i, j, k,bi,bj)-uVel(i, j, kp1,bi,bj)) +
673     $ (uVel(ip1,j, k,bi,bj)-uVel(ip1,j, kp1,bi,bj)) *
674     $ (uVel(ip1,j, k,bi,bj)-uVel(ip1,j, kp1,bi,bj)) +
675     $ (vVel(i, j, k,bi,bj)-vVel(i, j, kp1,bi,bj)) *
676     $ (vVel(i, j, k,bi,bj)-vVel(i, j, kp1,bi,bj)) +
677     $ (vVel(i, jp1,k,bi,bj)-vVel(i, jp1,kp1,bi,bj)) *
678     $ (vVel(i, jp1,k,bi,bj)-vVel(i, jp1,kp1,bi,bj)) )
679     #ifdef KPP_SMOOTH_SHSQ
680     shsq(i,j,k) = p5 * shsq(i,j,k) + p125 * (
681     $ (uVel(i, jm1,k,bi,bj)-uVel(i, jm1,kp1,bi,bj)) *
682     $ (uVel(i, jm1,k,bi,bj)-uVel(i, jm1,kp1,bi,bj)) +
683     $ (uVel(ip1,jm1,k,bi,bj)-uVel(ip1,jm1,kp1,bi,bj)) *
684     $ (uVel(ip1,jm1,k,bi,bj)-uVel(ip1,jm1,kp1,bi,bj)) +
685     $ (uVel(i, jp1,k,bi,bj)-uVel(i, jp1,kp1,bi,bj)) *
686     $ (uVel(i, jp1,k,bi,bj)-uVel(i, jp1,kp1,bi,bj)) +
687     $ (uVel(ip1,jp1,k,bi,bj)-uVel(ip1,jp1,kp1,bi,bj)) *
688     $ (uVel(ip1,jp1,k,bi,bj)-uVel(ip1,jp1,kp1,bi,bj)) +
689     $ (vVel(im1,j, k,bi,bj)-vVel(im1,j, kp1,bi,bj)) *
690     $ (vVel(im1,j, k,bi,bj)-vVel(im1,j, kp1,bi,bj)) +
691     $ (vVel(im1,jp1,k,bi,bj)-vVel(im1,jp1,kp1,bi,bj)) *
692     $ (vVel(im1,jp1,k,bi,bj)-vVel(im1,jp1,kp1,bi,bj)) +
693     $ (vVel(ip1,j, k,bi,bj)-vVel(ip1,j, kp1,bi,bj)) *
694     $ (vVel(ip1,j, k,bi,bj)-vVel(ip1,j, kp1,bi,bj)) +
695     $ (vVel(ip1,jp1,k,bi,bj)-vVel(ip1,jp1,kp1,bi,bj)) *
696     $ (vVel(ip1,jp1,k,bi,bj)-vVel(ip1,jp1,kp1,bi,bj)) )
697     #endif
698     END DO
699     END DO
700     END DO
701    
702     c shsq @ Nr computation
703     #ifdef FRUGAL_KPP
704     DO j = 1, sNy
705     DO i = 1, sNx
706     #else
707     DO j = jmin, jmax
708     DO i = imin, imax
709     #endif
710     shsq(i,j,Nr) = p0
711     END DO
712     END DO
713    
714     #ifndef FRUGAL_KPP
715     c set array edges to zero
716     DO k = 1, Nr
717     DO j = jmin, jmax
718     DO i = 1-OLx, imin-1
719     shsq(i,j,k) = p0
720     dVsq(i,j,k) = p0
721     END DO
722     DO i = imax+1, sNx+OLx
723     shsq(i,j,k) = p0
724     dVsq(i,j,k) = p0
725     END DO
726     END DO
727     DO i = 1-OLx, sNx+OLx
728     DO j = 1-OLy, jmin-1
729     shsq(i,j,k) = p0
730     dVsq(i,j,k) = p0
731     END DO
732     DO j = jmax+1, sNy+OLy
733     shsq(i,j,k) = p0
734     dVsq(i,j,k) = p0
735     END DO
736     END DO
737     END DO
738     #endif
739    
740     c-----------------------------------------------------------------------
741     c solve for viscosity, diffusivity, ghat, and hbl on "t-grid"
742     c-----------------------------------------------------------------------
743    
744     #ifdef FRUGAL_KPP
745     DO j = 1, sNy
746     DO i = 1, sNx
747     #else
748     DO j = 1-OLy, sNy+OLy
749     DO i = 1-OLx, sNx+OLx
750     #endif
751     work1(i,j) = nzmax(i,j,bi,bj)
752     work2(i,j) = Fcori(i,j,bi,bj)
753     END DO
754     END DO
755     CALL TIMER_START('KPPMIX [KPP_CALC]', myThid)
756     CALL KPPMIX (
757     I lri, work1, shsq, dVsq, ustar
758     I , bo, bosol, dbloc, Ritop, work2
759     I , ikey
760     O , vddiff
761     U , ghat
762     O , hbl
763     & )
764    
765     CALL TIMER_STOP ('KPPMIX [KPP_CALC]', myThid)
766    
767     IF (KPPmixingMaps) THEN
768     #ifdef FRUGAL_KPP
769     CALL PRINT_MAPRS(
770     I hbl, 'hbl', PRINT_MAP_XY,
771     I 1, sNx, 1, sNy, 1, 1, 1, 1,
772     I 1, sNx, 1, sNy, 1, -1, 1, 1, 1, 1, 1, 1, 1, 1, 1 )
773     #else
774     CALL PRINT_MAPRS(
775     I hbl, 'hbl', PRINT_MAP_XY,
776     I 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1, 1, 1, 1,
777     I 1, sNx, 1, sNy, 1, -1, 1, 1, 1, 1, 1, 1, 1, 1, 1 )
778     #endif
779     ENDIF
780    
781 heimbach 1.2 #ifdef ALLOW_AUTODIFF_TAMC
782     CADJ STORE vddiff, ghat = comlev1_kpp, key = ikey
783     #endif /* ALLOW_AUTODIFF_TAMC */
784 adcroft 1.1
785     c-----------------------------------------------------------------------
786     c zero out land values,
787     c make sure coefficients are within reasonable bounds,
788     c and transfer to global variables
789     c-----------------------------------------------------------------------
790    
791     #ifdef FRUGAL_KPP
792     DO j = 1, sNy
793     DO i = 1, sNx
794     #else
795     DO j = jmin, jmax
796     DO i = imin, imax
797     #endif
798     DO k = 1, Nr
799     c KPPviscAz
800     tempVar = min( maxKPPviscAz(k), vddiff(i,j,k-1,1) )
801     tempVar = max( minKPPviscAz, tempVar )
802     KPPviscAz(i,j,k,bi,bj) = tempVar*pMask(i,j,k,bi,bj)
803     c KPPdiffKzS
804     tempVar = min( maxKPPdiffKzS, vddiff(i,j,k-1,2) )
805     tempVar = max( minKPPdiffKzS, tempVar )
806     KPPdiffKzS(i,j,k,bi,bj) = tempVar*pMask(i,j,k,bi,bj)
807     c KPPdiffKzT
808     tempVar = min( maxKPPdiffKzT, vddiff(i,j,k-1,3) )
809     tempVar = max( minKPPdiffKzT, tempVar )
810     KPPdiffKzT(i,j,k,bi,bj) = tempVar*pMask(i,j,k,bi,bj)
811     c KPPghat
812     tempVar = min( maxKPPghat, ghat(i,j,k) )
813     tempVar = max( minKPPghat, tempVar )
814     KPPghat(i,j,k,bi,bj) = tempVar*pMask(i,j,k,bi,bj)
815     END DO
816     c KPPhbl: set to -zgrid(1) over land
817     KPPhbl(i,j,bi,bj) = (hbl(i,j) + zgrid(1))
818     & * pMask(i,j,1,bi,bj) -
819     & zgrid(1)
820     END DO
821     END DO
822     #ifdef FRUGAL_KPP
823     _EXCH_XYZ_R8(KPPviscAz , myThid )
824     _EXCH_XYZ_R8(KPPdiffKzS , myThid )
825     _EXCH_XYZ_R8(KPPdiffKzT , myThid )
826     _EXCH_XYZ_R8(KPPghat , myThid )
827     _EXCH_XY_R8 (KPPhbl , myThid )
828     #endif
829    
830     #ifdef KPP_SMOOTH_VISC
831     c horizontal smoothing of vertical viscosity
832     c as coded requires FRUGAL_KPP and OLx=4, OLy=4
833     c alternatively could recode with OLx=5, OLy=5
834    
835     DO k = 1, Nr
836     CALL SMOOTH_HORIZ_RL (
837     I k, bi, bj,
838     I KPPviscAz(1-OLx,1-OLy,k,bi,bj),
839     O KPPviscAz(1-OLx,1-OLy,k,bi,bj) )
840     END DO
841     #endif /* KPP_SMOOTH_VISC */
842    
843     #ifdef KPP_SMOOTH_DIFF
844     c horizontal smoothing of vertical diffusivity
845     c as coded requires FRUGAL_KPP and OLx=4, OLy=4
846     c alternatively could recode with OLx=5, OLy=5
847    
848     DO k = 1, Nr
849     CALL SMOOTH_HORIZ_RL (
850     I k, bi, bj,
851     I KPPdiffKzS(1-OLx,1-OLy,k,bi,bj),
852     O KPPdiffKzS(1-OLx,1-OLy,k,bi,bj) )
853     CALL SMOOTH_HORIZ_RL (
854     I k, bi, bj,
855     I KPPdiffKzT(1-OLx,1-OLy,k,bi,bj),
856     O KPPdiffKzT(1-OLx,1-OLy,k,bi,bj) )
857     END DO
858     #endif /* KPP_SMOOTH_DIFF */
859    
860    
861     C Compute fraction of solar short-wave flux penetrating to
862     C the bottom of the mixing layer.
863     DO j=1-OLy,sNy+OLy
864     DO i=1-OLx,sNx+OLx
865     worka(i,j) = KPPhbl(i,j,bi,bj)
866     ENDDO
867     ENDDO
868     CALL SWFRAC(
869 heimbach 1.2 I (sNx+2*OLx)*(sNy+2*OLy), m1, worka,
870 adcroft 1.1 O workb )
871     DO j=1-OLy,sNy+OLy
872     DO i=1-OLx,sNx+OLx
873     KPPfrac(i,j,bi,bj) = workb(i,j)
874     ENDDO
875     ENDDO
876    
877     ENDIF
878    
879     #endif ALLOW_KPP
880    
881     RETURN
882     END

  ViewVC Help
Powered by ViewVC 1.1.22