/[MITgcm]/MITgcm/pkg/ggl90/ggl90_calc.F
ViewVC logotype

Annotation of /MITgcm/pkg/ggl90/ggl90_calc.F

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


Revision 1.7 - (hide annotations) (download)
Tue Jun 6 22:15:19 2006 UTC (17 years, 11 months ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint58l_post, checkpoint58u_post, checkpoint58w_post, checkpoint60, checkpoint61, checkpoint58r_post, checkpoint58n_post, checkpoint58x_post, checkpoint58t_post, checkpoint58h_post, checkpoint58q_post, checkpoint59q, checkpoint59p, checkpoint59r, checkpoint58j_post, checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint59c, checkpoint59b, checkpoint59m, checkpoint59l, checkpoint59o, checkpoint59n, checkpoint59i, checkpoint59h, checkpoint59k, checkpoint59j, checkpoint59, checkpoint58i_post, checkpoint58g_post, checkpoint58o_post, checkpoint58y_post, checkpoint58k_post, checkpoint58v_post, checkpoint58s_post, checkpoint61b, checkpoint58p_post, checkpoint61a, checkpoint58m_post
Changes since 1.6: +8 -4 lines
 fix bug: avoid division by zero

1 mlosch 1.7 C $Header: /u/gcmpack/MITgcm/pkg/ggl90/ggl90_calc.F,v 1.6 2006/06/06 16:18:18 mlosch Exp $
2 mlosch 1.1 C $Name: $
3    
4     #include "GGL90_OPTIONS.h"
5    
6     CBOP
7     C !ROUTINE: GGL90_CALC
8    
9     C !INTERFACE: ======================================================
10     subroutine GGL90_CALC(
11     I bi, bj, myTime, myThid )
12    
13     C !DESCRIPTION: \bv
14     C /==========================================================\
15     C | SUBROUTINE GGL90_CALC |
16     C | o Compute all GGL90 fields defined in GGL90.h |
17     C |==========================================================|
18     C | Equation numbers refer to |
19     C | Gaspar et al. (1990), JGR 95 (C9), pp 16,179 |
20     C | Some parts of the implementation follow Blanke and |
21     C | Delecuse (1993), JPO, and OPA code, in particular the |
22     C | computation of the |
23     C | mixing length = max(min(lk,depth),lkmin) |
24     C \==========================================================/
25     IMPLICIT NONE
26     C
27     C--------------------------------------------------------------------
28    
29     C global parameters updated by ggl90_calc
30     C GGL90TKE - sub-grid turbulent kinetic energy (m^2/s^2)
31     C GGL90viscAz - GGL90 eddy viscosity coefficient (m^2/s)
32     C GGL90diffKzT - GGL90 diffusion coefficient for temperature (m^2/s)
33     C
34     C \ev
35    
36     C !USES: ============================================================
37     #include "SIZE.h"
38     #include "EEPARAMS.h"
39     #include "PARAMS.h"
40     #include "DYNVARS.h"
41     #include "GGL90.h"
42     #include "FFIELDS.h"
43     #include "GRID.h"
44    
45     C !INPUT PARAMETERS: ===================================================
46     c Routine arguments
47     c bi, bj - array indices on which to apply calculations
48     c myTime - Current time in simulation
49    
50     INTEGER bi, bj
51     INTEGER myThid
52     _RL myTime
53    
54     #ifdef ALLOW_GGL90
55    
56     C !LOCAL VARIABLES: ====================================================
57     c Local constants
58     C iMin, iMax, jMin, jMax, I, J - array computation indices
59     C K, Kp1, km1, kSurf, kBottom - vertical loop indices
60     C ab15, ab05 - weights for implicit timestepping
61     C uStarSquare - square of friction velocity
62     C verticalShear - (squared) vertical shear of horizontal velocity
63     C Nsquare - squared buoyancy freqency
64     C RiNumber - local Richardson number
65     C KappaM - (local) viscosity parameter (eq.10)
66     C KappaH - (local) diffusivity parameter for temperature (eq.11)
67     C KappaE - (local) diffusivity parameter for TKE (eq.15)
68     C buoyFreq - buoyancy freqency
69     C TKEdissipation - dissipation of TKE
70     C GGL90mixingLength- mixing length of scheme following Banke+Delecuse
71 mlosch 1.7 C rMixingLength- inverse of mixing length
72 mlosch 1.1 C totalDepth - thickness of water column (inverse of recip_Rcol)
73     C TKEPrandtlNumber - here, an empirical function of the Richardson number
74     C rhoK, rhoKm1 - density at layer K and Km1 (relative to K)
75     C gTKE - right hand side of implicit equation
76     INTEGER iMin ,iMax ,jMin ,jMax
77     INTEGER I, J, K, Kp1, Km1, kSurf, kBottom
78     _RL ab15, ab05
79     _RL uStarSquare
80     _RL verticalShear
81     _RL KappaM, KappaH
82     _RL Nsquare
83     _RL deltaTggl90
84     _RL SQRTTKE
85     _RL RiNumber
86     _RL TKEdissipation
87     _RL tempU, tempV, prTemp
88     _RL TKEPrandtlNumber (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
89     _RL GGL90mixingLength(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
90 mlosch 1.7 _RL rMixingLength(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
91 mlosch 1.1 _RL KappaE (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
92     _RL rhoK (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
93     _RL rhoKm1 (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94     _RL totalDepth (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95     _RL gTKE (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
96     C tri-diagonal matrix
97     _RL a(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
98     _RL b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
99     _RL c(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
100 mlosch 1.2 #ifdef ALLOW_GGL90_HORIZDIFF
101     C xA, yA - area of lateral faces
102     C dfx, dfy - diffusive flux across lateral faces
103     _RL xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
104     _RL yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
105     _RL dfx(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
106     _RL dfy(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
107     #endif /* ALLOW_GGL90_HORIZDIFF */
108 mlosch 1.1 CEOP
109     iMin = 2-OLx
110     iMax = sNx+OLx-1
111     jMin = 2-OLy
112     jMax = sNy+OLy-1
113    
114     C set separate time step (should be deltaTtracer)
115 jmc 1.4 deltaTggl90 = dTtracerLev(1)
116 mlosch 1.1 C
117     kSurf = 1
118     C implicit timestepping weights for dissipation
119     ab15 = 1.5 _d 0
120     ab05 = -0.5 _d 0
121     ab15 = 1. _d 0
122     ab05 = 0. _d 0
123    
124     C Initialize local fields
125     DO K = 1, Nr
126     DO J=1-Oly,sNy+Oly
127     DO I=1-Olx,sNx+Olx
128     gTKE(I,J,K) = 0. _d 0
129     KappaE(I,J,K) = 0. _d 0
130     TKEPrandtlNumber(I,J,K) = 0. _d 0
131 mlosch 1.7 GGL90mixingLength(I,J,K) = GGL90mixingLengthMin
132     rMixingLength(I,J,K) = 0. _d 0
133 mlosch 1.1 ENDDO
134     ENDDO
135     ENDDO
136     DO J=1-Oly,sNy+Oly
137     DO I=1-Olx,sNx+Olx
138     rhoK (I,J) = 0. _d 0
139     rhoKm1 (I,J) = 0. _d 0
140     totalDepth(I,J) = 0. _d 0
141     IF ( recip_Rcol(I,J,bi,bj) .NE. 0. )
142     & totalDepth(I,J) = 1./recip_Rcol(I,J,bi,bj)
143     ENDDO
144     ENDDO
145    
146     C start k-loop
147     DO K = 2, Nr
148     Km1 = K-1
149     Kp1 = MIN(Nr,K+1)
150     CALL FIND_RHO(
151     I bi, bj, iMin, iMax, jMin, jMax, Km1, K,
152     I theta, salt,
153     O rhoKm1,
154     I myThid )
155     CALL FIND_RHO(
156     I bi, bj, iMin, iMax, jMin, jMax, K, K,
157     I theta, salt,
158     O rhoK,
159     I myThid )
160     DO J=jMin,jMax
161     DO I=iMin,iMax
162     SQRTTKE=SQRT( GGL90TKE(I,J,K,bi,bj) )
163     C
164     C buoyancy frequency
165     C
166     Nsquare = - gravity*recip_rhoConst*recip_drC(K)
167     & * ( rhoKm1(I,J) - rhoK(I,J) )*maskC(I,J,K,bi,bj)
168     C vertical shear term (dU/dz)^2+(dV/dz)^2
169     tempu= .5*( uVel(I,J,Km1,bi,bj)+uVel(I+1,J,Km1,bi,bj)
170     & - (uVel(I,J,K ,bi,bj)+uVel(I+1,J,K ,bi,bj)) )
171     & *recip_drC(K)
172     tempv= .5*( vVel(I,J,Km1,bi,bj)+vVel(I,J+1,Km1,bi,bj)
173     & - (vVel(I,J,K ,bi,bj)+vVel(I,J+1,K ,bi,bj)) )
174     & *recip_drC(K)
175     verticalShear = tempU*tempU + tempV*tempV
176     RiNumber = MAX(Nsquare,0. _d 0)/(verticalShear+GGL90eps)
177     C compute Prandtl number (always greater than 0)
178     prTemp = 1. _d 0
179     IF ( RiNumber .GE. 0.2 ) prTemp = 5.0 * RiNumber
180 ce107 1.5 TKEPrandtlNumber(I,J,K) = MIN(10.0 _d 0,prTemp)
181 mlosch 1.1 C mixing length
182     GGL90mixingLength(I,J,K) =
183     & SQRTTKE/SQRT( MAX(Nsquare,GGL90eps) )
184     C impose upper bound for mixing length (total depth)
185     GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K),
186     & totalDepth(I,J))
187     C impose minimum mixing length (to avoid division by zero)
188     GGL90mixingLength(I,J,K) = MAX(GGL90mixingLength(I,J,K),
189     & GGL90mixingLengthMin)
190 mlosch 1.7 rMixingLength(I,J,K) = 1. _d 0 /GGL90mixingLength(I,J,K)
191 mlosch 1.1 C viscosity of last timestep
192     KappaM = GGL90ck*GGL90mixingLength(I,J,K)*SQRTTKE
193     KappaE(I,J,K) = KappaM*GGL90alpha
194     C dissipation term
195     TKEdissipation = ab05*GGL90ceps
196 mlosch 1.7 & *SQRTTKE*rMixingLength(I,J,K)
197 mlosch 1.1 & *GGL90TKE(I,J,K,bi,bj)
198     C sum up contributions to form the right hand side
199     gTKE(I,J,K) = GGL90TKE(I,J,K,bi,bj)
200     & + deltaTggl90*(
201     & + KappaM*verticalShear
202     & - KappaM*Nsquare/TKEPrandtlNumber(I,J,K)
203     & - TKEdissipation
204     & )
205     ENDDO
206     ENDDO
207 mlosch 1.2 ENDDO
208     C horizontal diffusion of TKE (requires an exchange in
209     C do_fields_blocking_exchanges)
210     #ifdef ALLOW_GGL90_HORIZDIFF
211     IF ( GGL90diffTKEh .GT. 0. _d 0 ) THEN
212     DO K=2,Nr
213     C common factors
214     DO j=1-Oly,sNy+Oly
215     DO i=1-Olx,sNx+Olx
216     xA(i,j) = _dyG(i,j,bi,bj)
217     & *drF(k)*_hFacW(i,j,k,bi,bj)
218     yA(i,j) = _dxG(i,j,bi,bj)
219     & *drF(k)*_hFacS(i,j,k,bi,bj)
220     ENDDO
221     ENDDO
222     C Compute diffusive fluxes
223     C ... across x-faces
224     DO j=1-Oly,sNy+Oly
225     dfx(1-Olx,j)=0.
226     DO i=1-Olx+1,sNx+Olx
227     dfx(i,j) = -GGL90diffTKEh*xA(i,j)
228     & *_recip_dxC(i,j,bi,bj)
229     & *(GGL90TKE(i,j,k,bi,bj)-GGL90TKE(i-1,j,k,bi,bj))
230     & *CosFacU(j,bi,bj)
231     ENDDO
232     ENDDO
233     C ... across y-faces
234     DO i=1-Olx,sNx+Olx
235     dfy(i,1-Oly)=0.
236     ENDDO
237     DO j=1-Oly+1,sNy+Oly
238     DO i=1-Olx,sNx+Olx
239     dfy(i,j) = -GGL90diffTKEh*yA(i,j)
240     & *_recip_dyC(i,j,bi,bj)
241     & *(GGL90TKE(i,j,k,bi,bj)-GGL90TKE(i,j-1,k,bi,bj))
242     #ifdef ISOTROPIC_COS_SCALING
243     & *CosFacV(j,bi,bj)
244     #endif /* ISOTROPIC_COS_SCALING */
245     ENDDO
246     ENDDO
247     C Compute divergence of fluxes
248     DO j=1-Oly,sNy+Oly-1
249     DO i=1-Olx,sNx+Olx-1
250     gTKE(i,j,k)=gTKE(i,j,k)
251     & -_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)*recip_rA(i,j,bi,bj)
252     & *( (dfx(i+1,j)-dfx(i,j))
253     & +(dfy(i,j+1)-dfy(i,j))
254     & )
255     ENDDO
256     ENDDO
257     C end of k-loop
258     ENDDO
259     C end if GGL90diffTKEh .eq. 0.
260     ENDIF
261     #endif /* ALLOW_GGL90_HORIZDIFF */
262    
263     C ============================================
264     C Implicit time step to update TKE for k=1,Nr;
265     C TKE(Nr+1)=0 by default
266     C ============================================
267 mlosch 1.1 C set up matrix
268 mlosch 1.2 C-- Lower diagonal
269 mlosch 1.1 DO j=jMin,jMax
270     DO i=iMin,iMax
271     a(i,j,1) = 0. _d 0
272     ENDDO
273     ENDDO
274     DO k=2,Nr
275     km1=MAX(1,k-1)
276     DO j=jMin,jMax
277     DO i=iMin,iMax
278     a(i,j,k) = -deltaTggl90
279     & *recip_drF(km1)*recip_hFacI(i,j,k,bi,bj)
280     & *.5*(KappaE(i,j, k )+KappaE(i,j,km1))
281     & *recip_drC(k)
282     IF (recip_hFacI(i,j,km1,bi,bj).EQ.0.) a(i,j,k)=0.
283     ENDDO
284     ENDDO
285     ENDDO
286 mlosch 1.2 C-- Upper diagonal
287 mlosch 1.1 DO j=jMin,jMax
288     DO i=iMin,iMax
289     c(i,j,1) = 0. _d 0
290     c(i,j,Nr) = 0. _d 0
291     ENDDO
292     ENDDO
293     CML DO k=1,Nr-1
294 mlosch 1.2 DO k=2,Nr-1
295 mlosch 1.1 kp1=min(Nr,k+1)
296     DO j=jMin,jMax
297     DO i=iMin,iMax
298     c(i,j,k) = -deltaTggl90
299     & *recip_drF( k )*recip_hFacI(i,j,k,bi,bj)
300     & *.5*(KappaE(i,j,k)+KappaE(i,j,kp1))
301     & *recip_drC(k)
302 mlosch 1.2 IF (recip_hFacI(i,j,kp1,bi,bj).EQ.0.) c(i,j,k)=0.
303 mlosch 1.1 ENDDO
304     ENDDO
305     ENDDO
306 mlosch 1.2 C-- Center diagonal
307 mlosch 1.1 DO k=1,Nr
308     DO j=jMin,jMax
309     DO i=iMin,iMax
310     b(i,j,k) = 1. _d 0 - c(i,j,k) - a(i,j,k)
311     & + ab15*deltaTggl90*GGL90ceps*SQRT(GGL90TKE(I,J,K,bi,bj))
312 mlosch 1.7 & *rMixingLength(I,J,K)
313 mlosch 1.1 ENDDO
314     ENDDO
315     ENDDO
316     C end set up matrix
317    
318     C
319     C Apply boundary condition
320     C
321     DO J=jMin,jMax
322     DO I=iMin,iMax
323     C estimate friction velocity uStar from surface forcing
324     uStarSquare = SQRT(
325     & ( .5*( surfaceForcingU(I, J, bi,bj)
326     & + surfaceForcingU(I+1,J, bi,bj) ) )**2
327     & + ( .5*( surfaceForcingV(I, J, bi,bj)
328     & + surfaceForcingV(I, J+1,bi,bj) ) )**2
329     & )
330     C Dirichlet surface boundary condition for TKE
331 mlosch 1.6 gTKE(I,J,kSurf) = MAX(GGL90TKEsurfMin,GGL90m2*uStarSquare)
332 mlosch 1.1 & *maskC(I,J,kSurf,bi,bj)
333     C Dirichlet bottom boundary condition for TKE = GGL90TKEbottom
334     kBottom = MIN(MAX(kLowC(I,J,bi,bj),1),Nr)
335     gTKE(I,J,kBottom) = gTKE(I,J,kBottom)
336     & - GGL90TKEbottom*c(I,J,kBottom)
337     ENDDO
338     ENDDO
339     C
340     C solve tri-diagonal system, and store solution on gTKE (previously rhs)
341     C
342     CALL GGL90_SOLVE( bi, bj, iMin, iMax, jMin, jMax,
343     I a, b, c,
344     U gTKE,
345     I myThid )
346     C
347     C now update TKE
348     C
349     DO K=1,Nr
350     DO J=jMin,jMax
351     DO I=iMin,iMax
352     C impose minimum TKE to avoid numerical undershoots below zero
353     GGL90TKE(I,J,K,bi,bj) = MAX( gTKE(I,J,K), GGL90TKEmin )
354     & * maskC(I,J,K,bi,bj)
355     ENDDO
356     ENDDO
357     ENDDO
358     C
359 mlosch 1.2 C end of time step
360     C ===============================
361 mlosch 1.1 C compute viscosity coefficients
362     C
363     DO K=2,Nr
364     DO J=jMin,jMax
365     DO I=iMin,iMax
366     C Eq. (11), (18) and (21)
367     KappaM = GGL90ck*GGL90mixingLength(I,J,K)*
368     & SQRT( GGL90TKE(I,J,K,bi,bj) )
369     KappaH = KappaM/TKEPrandtlNumber(I,J,K)
370     C Set a minium (= background) value
371     KappaM = MAX(KappaM,viscAr)
372 jmc 1.3 KappaH = MAX(KappaH,diffKrNrT(k))
373 mlosch 1.1 C Set a maximum and mask land point
374     GGL90viscAr(I,J,K,bi,bj) = MIN(KappaM,GGL90viscMax)
375     & * maskC(I,J,K,bi,bj)
376     GGL90diffKr(I,J,K,bi,bj) = MIN(KappaH,GGL90diffMax)
377     & * maskC(I,J,K,bi,bj)
378     ENDDO
379     ENDDO
380     C end third k-loop
381     ENDDO
382    
383     #endif /* ALLOW_GGL90 */
384    
385     RETURN
386     END
387    

  ViewVC Help
Powered by ViewVC 1.1.22