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

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

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

revision 1.3 by jmc, Sun Oct 17 23:11:59 2004 UTC revision 1.16 by gforget, Fri Aug 6 18:37:05 2010 UTC
# Line 11  C !INTERFACE: ========================== Line 11  C !INTERFACE: ==========================
11       I     bi, bj, myTime, myThid )       I     bi, bj, myTime, myThid )
12    
13  C !DESCRIPTION: \bv  C !DESCRIPTION: \bv
14  C     /==========================================================\  C     *==========================================================*
15  C     | SUBROUTINE GGL90_CALC                                    |  C     | SUBROUTINE GGL90_CALC                                    |
16  C     | o Compute all GGL90 fields defined in GGL90.h            |  C     | o Compute all GGL90 fields defined in GGL90.h            |
17  C     |==========================================================|  C     *==========================================================*
18  C     | Equation numbers refer to                                |  C     | Equation numbers refer to                                |
19  C     | Gaspar et al. (1990), JGR 95 (C9), pp 16,179             |  C     | Gaspar et al. (1990), JGR 95 (C9), pp 16,179             |
20  C     | Some parts of the implementation follow Blanke and       |  C     | Some parts of the implementation follow Blanke and       |
21  C     | Delecuse (1993), JPO, and OPA code, in particular the    |  C     | Delecuse (1993), JPO, and OPA code, in particular the    |
22  C     | computation of the                                       |  C     | computation of the                                       |
23  C     | mixing length = max(min(lk,depth),lkmin)                 |  C     | mixing length = max(min(lk,depth),lkmin)                 |
24  C     \==========================================================/  C     *==========================================================*
       IMPLICIT NONE  
 C  
 C--------------------------------------------------------------------  
25    
26  C global parameters updated by ggl90_calc  C global parameters updated by ggl90_calc
27  C     GGL90TKE     - sub-grid turbulent kinetic energy           (m^2/s^2)  C     GGL90TKE     :: sub-grid turbulent kinetic energy          (m^2/s^2)
28  C     GGL90viscAz  - GGL90 eddy viscosity coefficient              (m^2/s)  C     GGL90viscAz  :: GGL90 eddy viscosity coefficient             (m^2/s)
29  C     GGL90diffKzT - GGL90 diffusion coefficient for temperature   (m^2/s)  C     GGL90diffKzT :: GGL90 diffusion coefficient for temperature  (m^2/s)
 C  
30  C \ev  C \ev
31    
32  C !USES: ============================================================  C !USES: ============================================================
33          IMPLICIT NONE
34  #include "SIZE.h"  #include "SIZE.h"
35  #include "EEPARAMS.h"  #include "EEPARAMS.h"
36  #include "PARAMS.h"  #include "PARAMS.h"
# Line 43  C !USES: =============================== Line 40  C !USES: ===============================
40  #include "GRID.h"  #include "GRID.h"
41    
42  C !INPUT PARAMETERS: ===================================================  C !INPUT PARAMETERS: ===================================================
43  c Routine arguments  C Routine arguments
44  c     bi, bj - array indices on which to apply calculations  C     bi, bj :: array indices on which to apply calculations
45  c     myTime - Current time in simulation  C     myTime :: Current time in simulation
46    C     myThid :: My Thread Id number
47        INTEGER bi, bj        INTEGER bi, bj
       INTEGER myThid  
48        _RL     myTime        _RL     myTime
49          INTEGER myThid
50    CEOP
51    
52  #ifdef ALLOW_GGL90  #ifdef ALLOW_GGL90
53    
54  C !LOCAL VARIABLES: ====================================================  C !LOCAL VARIABLES: ====================================================
55  c Local constants  C Local constants
56  C     iMin, iMax, jMin, jMax, I, J - array computation indices  C     iMin, iMax, jMin, jMax, I, J - array computation indices
57  C     K, Kp1, km1, kSurf, kBottom  - vertical loop indices  C     K, Kp1, km1, kSurf, kBottom  - vertical loop indices
58  C     ab15, ab05       - weights for implicit timestepping  C     ab15, ab05       - weights for implicit timestepping
59  C     uStarSquare      - square of friction velocity  C     uStarSquare      - square of friction velocity
60  C     verticalShear    - (squared) vertical shear of horizontal velocity  C     verticalShear    - (squared) vertical shear of horizontal velocity
61  C     Nsquare          - squared buoyancy freqency  C     Nsquare          - squared buoyancy freqency
62  C     RiNumber         - local Richardson number  C     RiNumber         - local Richardson number
63  C     KappaM           - (local) viscosity parameter (eq.10)  C     KappaM           - (local) viscosity parameter (eq.10)
64  C     KappaH           - (local) diffusivity parameter for temperature (eq.11)  C     KappaH           - (local) diffusivity parameter for temperature (eq.11)
65  C     KappaE           - (local) diffusivity parameter for TKE (eq.15)  C     KappaE           - (local) diffusivity parameter for TKE (eq.15)
 C     buoyFreq         - buoyancy freqency  
66  C     TKEdissipation   - dissipation of TKE  C     TKEdissipation   - dissipation of TKE
67  C     GGL90mixingLength- mixing length of scheme following Banke+Delecuse  C     GGL90mixingLength- mixing length of scheme following Banke+Delecuse
68    C         rMixingLength- inverse of mixing length
69  C     totalDepth       - thickness of water column (inverse of recip_Rcol)  C     totalDepth       - thickness of water column (inverse of recip_Rcol)
70  C     TKEPrandtlNumber - here, an empirical function of the Richardson number  C     TKEPrandtlNumber - here, an empirical function of the Richardson number
71  C     rhoK, rhoKm1     - density at layer K and Km1 (relative to K)  C     rhoK, rhoKm1     - density at layer K and Km1 (relative to K)
# Line 78  C     gTKE             - right hand side Line 76  C     gTKE             - right hand side
76        _RL     uStarSquare        _RL     uStarSquare
77        _RL     verticalShear        _RL     verticalShear
78        _RL     KappaM, KappaH        _RL     KappaM, KappaH
79        _RL     Nsquare  c     _RL     Nsquare
80          _RL     Nsquare(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
81        _RL     deltaTggl90        _RL     deltaTggl90
82        _RL     SQRTTKE  c     _RL     SQRTTKE
83          _RL     SQRTTKE(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
84        _RL     RiNumber        _RL     RiNumber
85        _RL     TKEdissipation        _RL     TKEdissipation
86        _RL     tempU, tempV, prTemp        _RL     tempU, tempV, prTemp
87          _RL     MaxLength, tmpmlx, tmpVisc
88        _RL     TKEPrandtlNumber (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL     TKEPrandtlNumber (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
89        _RL     GGL90mixingLength(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL     GGL90mixingLength(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
90          _RL     rMixingLength    (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
91          _RL     mxLength_Dn      (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
92        _RL     KappaE           (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL     KappaE           (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
93        _RL     rhoK             (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL     rhoK             (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
94        _RL     rhoKm1           (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL     rhoKm1           (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
95        _RL     totalDepth       (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL     totalDepth       (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
96        _RL     gTKE             (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL     gTKE             (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
97  C     tri-diagonal matrix        _RL     GGL90visctmp     (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
98    C-    tri-diagonal matrix
99        _RL     a(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL     a(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
100        _RL     b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL     b(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
101        _RL     c(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)        _RL     c(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
102          INTEGER errCode
103  #ifdef ALLOW_GGL90_HORIZDIFF  #ifdef ALLOW_GGL90_HORIZDIFF
104  C     xA, yA   - area of lateral faces  C-    xA, yA   - area of lateral faces
105  C     dfx, dfy - diffusive flux across lateral faces  C-    dfx, dfy - diffusive flux across lateral faces
106        _RL     xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL     xA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
107        _RL     yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL     yA (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
108        _RL     dfx(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL     dfx(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
109        _RL     dfy(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL     dfy(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
110  #endif /* ALLOW_GGL90_HORIZDIFF */  #endif /* ALLOW_GGL90_HORIZDIFF */
111  CEOP  #ifdef ALLOW_GGL90_SMOOTH
112          _RL p4, p8, p16
113          p4=0.25 _d 0
114          p8=0.125 _d 0
115          p16=0.0625 _d 0
116    #endif
117        iMin = 2-OLx        iMin = 2-OLx
118        iMax = sNx+OLx-1        iMax = sNx+OLx-1
119        jMin = 2-OLy        jMin = 2-OLy
120        jMax = sNy+OLy-1        jMax = sNy+OLy-1
121    
122  C     set separate time step (should be deltaTtracer)  C     set separate time step (should be deltaTtracer)
123        deltaTggl90 = deltaTtracer        deltaTggl90 = dTtracerLev(1)
124  C      
125        kSurf = 1        kSurf = 1
126  C     implicit timestepping weights for dissipation  C     implicit timestepping weights for dissipation
127        ab15 =  1.5 _d 0        ab15 =  1.5 _d 0
# Line 125  C     Initialize local fields Line 135  C     Initialize local fields
135          DO I=1-Olx,sNx+Olx          DO I=1-Olx,sNx+Olx
136           gTKE(I,J,K)              = 0. _d 0           gTKE(I,J,K)              = 0. _d 0
137           KappaE(I,J,K)            = 0. _d 0           KappaE(I,J,K)            = 0. _d 0
138           TKEPrandtlNumber(I,J,K)  = 0. _d 0           TKEPrandtlNumber(I,J,K)  = 1. _d 0
139           GGL90mixingLength(I,J,K) = 0. _d 0           GGL90mixingLength(I,J,K) = GGL90mixingLengthMin
140          ENDDO          ENDDO
141         ENDDO             ENDDO
142        ENDDO        ENDDO
143        DO J=1-Oly,sNy+Oly        DO J=1-Oly,sNy+Oly
144         DO I=1-Olx,sNx+Olx         DO I=1-Olx,sNx+Olx
145          rhoK    (I,J) = 0. _d 0          rhoK(I,J)          = 0. _d 0
146          rhoKm1  (I,J) = 0. _d 0          rhoKm1(I,J)        = 0. _d 0
147          totalDepth(I,J)   = 0. _d 0          totalDepth(I,J)    = Ro_surf(i,j,bi,bj) - R_low(i,j,bi,bj)
148          IF ( recip_Rcol(I,J,bi,bj) .NE. 0. )          rMixingLength(i,j,1) = 0. _d 0
149       &       totalDepth(I,J) = 1./recip_Rcol(I,J,bi,bj)          mxLength_Dn(I,J,1) = GGL90mixingLengthMin
150            SQRTTKE(i,j,1) = SQRT( GGL90TKE(i,j,1,bi,bj) )
151         ENDDO         ENDDO
152        ENDDO        ENDDO
153    
154  C     start k-loop  C     start k-loop
155        DO K = 2, Nr        DO K = 2, Nr
156         Km1 = K-1         Km1 = K-1
157         Kp1 = MIN(Nr,K+1)  c      Kp1 = MIN(Nr,K+1)
158         CALL FIND_RHO(         CALL FIND_RHO_2D(
159       I      bi, bj, iMin, iMax, jMin, jMax, Km1, K,       I      iMin, iMax, jMin, jMax, K,
160       I      theta, salt,       I      theta(1-OLx,1-OLy,Km1,bi,bj), salt(1-OLx,1-OLy,Km1,bi,bj),
161       O      rhoKm1,       O      rhoKm1,
162       I      myThid )       I      Km1, bi, bj, myThid )
163         CALL FIND_RHO(  
164       I      bi, bj, iMin, iMax, jMin, jMax, K, K,         CALL FIND_RHO_2D(
165       I      theta, salt,       I      iMin, iMax, jMin, jMax, K,
166         I      theta(1-OLx,1-OLy,K,bi,bj), salt(1-OLx,1-OLy,K,bi,bj),
167       O      rhoK,       O      rhoK,
168       I      myThid )       I      K, bi, bj, myThid )
169         DO J=jMin,jMax         DO J=jMin,jMax
170          DO I=iMin,iMax          DO I=iMin,iMax
171           SQRTTKE=SQRT( GGL90TKE(I,J,K,bi,bj) )           SQRTTKE(i,j,k)=SQRT( GGL90TKE(I,J,K,bi,bj) )
172  C  
173  C     buoyancy frequency  C     buoyancy frequency
174  C           Nsquare(i,j,k) = - gravity*recip_rhoConst*recip_drC(K)
          Nsquare = - gravity*recip_rhoConst*recip_drC(K)  
175       &        * ( rhoKm1(I,J) - rhoK(I,J) )*maskC(I,J,K,bi,bj)       &        * ( rhoKm1(I,J) - rhoK(I,J) )*maskC(I,J,K,bi,bj)
176    cC     vertical shear term (dU/dz)^2+(dV/dz)^2
177    c         tempU= .5 _d 0*( uVel(I,J,Km1,bi,bj)+uVel(I+1,J,Km1,bi,bj)
178    c     &                 -( uVel(I,J,K  ,bi,bj)+uVel(I+1,J,K  ,bi,bj)) )
179    c     &        *recip_drC(K)
180    c         tempV= .5 _d 0*( vVel(I,J,Km1,bi,bj)+vVel(I,J+1,Km1,bi,bj)
181    c     &                 -( vVel(I,J,K  ,bi,bj)+vVel(I,J+1,K  ,bi,bj)) )
182    c     &        *recip_drC(K)
183    c         verticalShear = tempU*tempU + tempV*tempV
184    c         RiNumber   = MAX(Nsquare(i,j,k),0. _d 0)/(verticalShear+GGL90eps)
185    cC     compute Prandtl number (always greater than 0)
186    c         prTemp = 1. _d 0
187    c         IF ( RiNumber .GE. 0.2 _d 0 ) prTemp = 5. _d 0 * RiNumber
188    c         TKEPrandtlNumber(I,J,K) = MIN(10. _d 0,prTemp)
189    C     mixing length
190             GGL90mixingLength(I,J,K) = SQRTTWO *
191         &        SQRTTKE(i,j,k)/SQRT( MAX(Nsquare(i,j,k),GGL90eps) )
192            ENDDO
193           ENDDO
194          ENDDO
195    
196    C- Impose upper and lower bound for mixing length
197          IF ( mxlMaxFlag .EQ. 0 ) THEN
198    C-
199           DO k=2,Nr
200            DO J=jMin,jMax
201             DO I=iMin,iMax
202              MaxLength=totalDepth(I,J)
203              GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K),
204         &                                   MaxLength)
205             ENDDO
206            ENDDO
207           ENDDO
208    
209           DO k=2,Nr
210            DO J=jMin,jMax
211             DO I=iMin,iMax
212              GGL90mixingLength(I,J,K) = MAX(GGL90mixingLength(I,J,K),
213         &                                   GGL90mixingLengthMin)
214              rMixingLength(I,J,K) = 1. _d 0 / GGL90mixingLength(I,J,K)
215             ENDDO
216            ENDDO
217           ENDDO
218    
219          ELSEIF ( mxlMaxFlag .EQ. 1 ) THEN
220    C-
221           DO k=2,Nr
222            DO J=jMin,jMax
223             DO I=iMin,iMax
224              MaxLength=MIN(Ro_surf(I,J,bi,bj)-rF(k),rF(k)-R_low(I,J,bi,bj))
225    c         MaxLength=MAX(MaxLength,20. _d 0)
226              GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K),
227         &                                   MaxLength)
228             ENDDO
229            ENDDO
230           ENDDO
231    
232           DO k=2,Nr
233            DO J=jMin,jMax
234             DO I=iMin,iMax
235              GGL90mixingLength(I,J,K) = MAX(GGL90mixingLength(I,J,K),
236         &                                   GGL90mixingLengthMin)
237              rMixingLength(I,J,K) = 1. _d 0 / GGL90mixingLength(I,J,K)
238             ENDDO
239            ENDDO
240           ENDDO
241    
242          ELSEIF ( mxlMaxFlag .EQ. 2 ) THEN
243    C-
244    cgf ensure mixing between first and second level
245    c      DO J=jMin,jMax
246    c        DO I=iMin,iMax
247    c         GGL90mixingLength(I,J,2)=drF(1)
248    c        ENDDO
249    c      ENDDO
250    cgf
251           DO k=2,Nr
252            DO J=jMin,jMax
253             DO I=iMin,iMax
254              GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K),
255         &        GGL90mixingLength(I,J,K-1)+drF(k-1))
256             ENDDO
257            ENDDO
258           ENDDO
259           DO J=jMin,jMax
260            DO I=iMin,iMax
261              GGL90mixingLength(I,J,Nr) = MIN(GGL90mixingLength(I,J,Nr),
262         &       GGL90mixingLengthMin+drF(Nr))
263            ENDDO
264           ENDDO
265           DO k=Nr-1,2,-1
266            DO J=jMin,jMax
267             DO I=iMin,iMax
268              GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K),
269         &        GGL90mixingLength(I,J,K+1)+drF(k))
270             ENDDO
271            ENDDO
272           ENDDO
273    
274           DO k=2,Nr
275            DO J=jMin,jMax
276             DO I=iMin,iMax
277              GGL90mixingLength(I,J,K) = MAX(GGL90mixingLength(I,J,K),
278         &                                   GGL90mixingLengthMin)
279              rMixingLength(I,J,K) = 1. _d 0 / GGL90mixingLength(I,J,K)
280             ENDDO
281            ENDDO
282           ENDDO
283    
284          ELSEIF ( mxlMaxFlag .EQ. 3 ) THEN
285    C-
286           DO k=2,Nr
287            DO J=jMin,jMax
288             DO I=iMin,iMax
289              mxLength_Dn(I,J,K) = MIN(GGL90mixingLength(I,J,K),
290         &        mxLength_Dn(I,J,K-1)+drF(k-1))
291             ENDDO
292            ENDDO
293           ENDDO
294           DO J=jMin,jMax
295            DO I=iMin,iMax
296              GGL90mixingLength(I,J,Nr) = MIN(GGL90mixingLength(I,J,Nr),
297         &       GGL90mixingLengthMin+drF(Nr))
298            ENDDO
299           ENDDO
300           DO k=Nr-1,2,-1
301            DO J=jMin,jMax
302             DO I=iMin,iMax
303              GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K),
304         &        GGL90mixingLength(I,J,K+1)+drF(k))
305             ENDDO
306            ENDDO
307           ENDDO
308    
309           DO k=2,Nr
310            DO J=jMin,jMax
311             DO I=iMin,iMax
312              GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K),
313         &                                  mxLength_Dn(I,J,K))
314              tmpmlx = SQRT( GGL90mixingLength(I,J,K)*mxLength_Dn(I,J,K) )
315              tmpmlx = MAX( tmpmlx, GGL90mixingLengthMin)
316              rMixingLength(I,J,K) = 1. _d 0 / tmpmlx
317             ENDDO
318            ENDDO
319           ENDDO
320    
321          ELSE
322           STOP 'GGL90_CALC: Wrong mxlMaxFlag (mixing length limit)'
323          ENDIF
324    
325    C- Impose minimum mixing length (to avoid division by zero)
326    c      DO k=2,Nr
327    c      DO J=jMin,jMax
328    c       DO I=iMin,iMax
329    c        GGL90mixingLength(I,J,K) = MAX(GGL90mixingLength(I,J,K),
330    c    &        GGL90mixingLengthMin)
331    c        rMixingLength(I,J,K) = 1. _d 0 /GGL90mixingLength(I,J,K)
332    c       ENDDO
333    c      ENDDO
334    c     ENDDO
335    
336    
337          DO k=2,Nr
338           Km1 = K-1
339           DO J=jMin,jMax
340            DO I=iMin,iMax
341  C     vertical shear term (dU/dz)^2+(dV/dz)^2  C     vertical shear term (dU/dz)^2+(dV/dz)^2
342           tempu= .5*(  uVel(I,J,Km1,bi,bj)+uVel(I+1,J,Km1,bi,bj)           tempU= .5 _d 0*( uVel(I,J,Km1,bi,bj)+uVel(I+1,J,Km1,bi,bj)
343       &             - (uVel(I,J,K  ,bi,bj)+uVel(I+1,J,K  ,bi,bj)) )       &                 -( uVel(I,J,K  ,bi,bj)+uVel(I+1,J,K  ,bi,bj)) )
344       &        *recip_drC(K)       &        *recip_drC(K)
345           tempv= .5*(  vVel(I,J,Km1,bi,bj)+vVel(I,J+1,Km1,bi,bj)           tempV= .5 _d 0*( vVel(I,J,Km1,bi,bj)+vVel(I,J+1,Km1,bi,bj)
346       &             - (vVel(I,J,K  ,bi,bj)+vVel(I,J+1,K  ,bi,bj)) )       &                 -( vVel(I,J,K  ,bi,bj)+vVel(I,J+1,K  ,bi,bj)) )
347       &        *recip_drC(K)       &        *recip_drC(K)
348           verticalShear = tempU*tempU + tempV*tempV           verticalShear = tempU*tempU + tempV*tempV
349           RiNumber   = MAX(Nsquare,0. _d 0)/(verticalShear+GGL90eps)           RiNumber = MAX(Nsquare(i,j,k),0. _d 0)/(verticalShear+GGL90eps)
350  C     compute Prandtl number (always greater than 0)  C     compute Prandtl number (always greater than 0)
351           prTemp = 1. _d 0           prTemp = 1. _d 0
352           IF ( RiNumber .GE. 0.2 ) prTemp = 5.0 * RiNumber           IF ( RiNumber .GE. 0.2 _d 0 ) prTemp = 5. _d 0 * RiNumber
353           TKEPrandtlNumber(I,J,K) = MIN(10.,prTemp)           TKEPrandtlNumber(I,J,K) = MIN(10. _d 0,prTemp)
354  C     mixing length  c         TKEPrandtlNumber(I,J,K) = 1. _d 0
355           GGL90mixingLength(I,J,K) =  
356       &        SQRTTKE/SQRT( MAX(Nsquare,GGL90eps) )  C     viscosity and diffusivity
357  C     impose upper bound for mixing length (total depth)           KappaM = GGL90ck*GGL90mixingLength(I,J,K)*SQRTTKE(i,j,k)
358           GGL90mixingLength(I,J,K) = MIN(GGL90mixingLength(I,J,K),           GGL90visctmp(I,J,K) = MAX(KappaM,diffKrNrT(k))
359       &        totalDepth(I,J))       &                            * maskC(I,J,K,bi,bj)
360  C     impose minimum mixing length (to avoid division by zero)  c        note: storing GGL90visctmp like this, and using it later to compute
361           GGL90mixingLength(I,J,K) = MAX(GGL90mixingLength(I,J,K),  c              GGL9rdiffKr etc. is robust in case of smoothing (e.g. see OPA)
362       &        GGL90mixingLengthMin)           KappaM = MAX(KappaM,viscArNr(k)) * maskC(I,J,K,bi,bj)
363  C     viscosity of last timestep           KappaH = KappaM/TKEPrandtlNumber(I,J,K)
364           KappaM = GGL90ck*GGL90mixingLength(I,J,K)*SQRTTKE           KappaE(I,J,K) = GGL90alpha * KappaM * maskC(I,J,K,bi,bj)
365           KappaE(I,J,K) = KappaM*GGL90alpha  
366  C     dissipation term  C     dissipation term
367           TKEdissipation = ab05*GGL90ceps           TKEdissipation = ab05*GGL90ceps
368       &        *SQRTTKE/GGL90mixingLength(I,J,K)       &        *SQRTTKE(i,j,k)*rMixingLength(I,J,K)
369       &        *GGL90TKE(I,J,K,bi,bj)             &        *GGL90TKE(I,J,K,bi,bj)
370  C     sum up contributions to form the right hand side  C     sum up contributions to form the right hand side
371           gTKE(I,J,K) = GGL90TKE(I,J,K,bi,bj)           gTKE(I,J,K) = GGL90TKE(I,J,K,bi,bj)
372       &        + deltaTggl90*(       &        + deltaTggl90*(
373       &        + KappaM*verticalShear       &        + KappaM*verticalShear
374       &        - KappaM*Nsquare/TKEPrandtlNumber(I,J,K)       &        - KappaH*Nsquare(i,j,k)
375       &        - TKEdissipation       &        - TKEdissipation
376       &        )       &        )
377          ENDDO            ENDDO
378         ENDDO         ENDDO
379        ENDDO        ENDDO
380  C     horizontal diffusion of TKE (requires an exchange in  
381  C       do_fields_blocking_exchanges)  C     horizontal diffusion of TKE (requires an exchange in
382    C      do_fields_blocking_exchanges)
383  #ifdef ALLOW_GGL90_HORIZDIFF  #ifdef ALLOW_GGL90_HORIZDIFF
384        IF ( GGL90diffTKEh .GT. 0. _d 0 ) THEN        IF ( GGL90diffTKEh .GT. 0. _d 0 ) THEN
385         DO K=2,Nr         DO K=2,Nr
# Line 214  C     common factors Line 391  C     common factors
391            yA(i,j) = _dxG(i,j,bi,bj)            yA(i,j) = _dxG(i,j,bi,bj)
392       &         *drF(k)*_hFacS(i,j,k,bi,bj)       &         *drF(k)*_hFacS(i,j,k,bi,bj)
393           ENDDO           ENDDO
394          ENDDO            ENDDO
395  C     Compute diffusive fluxes  C     Compute diffusive fluxes
396  C     ... across x-faces  C     ... across x-faces
397          DO j=1-Oly,sNy+Oly          DO j=1-Oly,sNy+Oly
398           dfx(1-Olx,j)=0.           dfx(1-Olx,j)=0. _d 0
399           DO i=1-Olx+1,sNx+Olx           DO i=1-Olx+1,sNx+Olx
400            dfx(i,j) = -GGL90diffTKEh*xA(i,j)            dfx(i,j) = -GGL90diffTKEh*xA(i,j)
401       &      *_recip_dxC(i,j,bi,bj)       &      *_recip_dxC(i,j,bi,bj)
# Line 228  C     ... across x-faces Line 405  C     ... across x-faces
405          ENDDO          ENDDO
406  C     ... across y-faces  C     ... across y-faces
407          DO i=1-Olx,sNx+Olx          DO i=1-Olx,sNx+Olx
408           dfy(i,1-Oly)=0.           dfy(i,1-Oly)=0. _d 0
409          ENDDO          ENDDO
410          DO j=1-Oly+1,sNy+Oly          DO j=1-Oly+1,sNy+Oly
411           DO i=1-Olx,sNx+Olx           DO i=1-Olx,sNx+Olx
# Line 239  C     ... across y-faces Line 416  C     ... across y-faces
416       &      *CosFacV(j,bi,bj)       &      *CosFacV(j,bi,bj)
417  #endif /* ISOTROPIC_COS_SCALING */  #endif /* ISOTROPIC_COS_SCALING */
418           ENDDO           ENDDO
419          ENDDO            ENDDO
420  C     Compute divergence of fluxes  C     Compute divergence of fluxes
421          DO j=1-Oly,sNy+Oly-1          DO j=1-Oly,sNy+Oly-1
422           DO i=1-Olx,sNx+Olx-1           DO i=1-Olx,sNx+Olx-1
423            gTKE(i,j,k)=gTKE(i,j,k)            gTKE(i,j,k)=gTKE(i,j,k)
424       &   -_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)*recip_rA(i,j,bi,bj)       &   -_recip_hFacC(i,j,k,bi,bj)*recip_drF(k)*recip_rA(i,j,bi,bj)
425       &         *( (dfx(i+1,j)-dfx(i,j))       &         *( (dfx(i+1,j)-dfx(i,j))
426       &           +(dfy(i,j+1)-dfy(i,j))       &           +(dfy(i,j+1)-dfy(i,j))
427       &           )       &           )*deltaTggl90
428           ENDDO             ENDDO
429          ENDDO          ENDDO
430  C       end of k-loop  C       end of k-loop
431         ENDDO         ENDDO
432  C     end if GGL90diffTKEh .eq. 0.  C     end if GGL90diffTKEh .eq. 0.
433        ENDIF        ENDIF
# Line 264  C     set up matrix Line 441  C     set up matrix
441  C--   Lower diagonal  C--   Lower diagonal
442        DO j=jMin,jMax        DO j=jMin,jMax
443         DO i=iMin,iMax         DO i=iMin,iMax
444           a(i,j,1) = 0. _d 0           a(i,j,1) = 0. _d 0
445         ENDDO         ENDDO
446        ENDDO        ENDDO
447        DO k=2,Nr        DO k=2,Nr
448         km1=MAX(1,k-1)         km1=MAX(2,k-1)
449         DO j=jMin,jMax         DO j=jMin,jMax
450          DO i=iMin,iMax          DO i=iMin,iMax
451    C-    We keep recip_hFacC in the diffusive flux calculation,
452    C-    but no hFacC in TKE volume control
453    C-    No need for maskC(k-1) with recip_hFacC(k-1)
454           a(i,j,k) = -deltaTggl90           a(i,j,k) = -deltaTggl90
455       &        *recip_drF(km1)*recip_hFacI(i,j,k,bi,bj)       &        *recip_drF(k-1)*recip_hFacC(i,j,k-1,bi,bj)
456       &        *.5*(KappaE(i,j, k )+KappaE(i,j,km1))       &        *.5 _d 0*(KappaE(i,j, k )+KappaE(i,j,km1))
457       &        *recip_drC(k)       &        *recip_drC(k)*maskC(i,j,k,bi,bj)
           IF (recip_hFacI(i,j,km1,bi,bj).EQ.0.) a(i,j,k)=0.  
458          ENDDO          ENDDO
459         ENDDO         ENDDO
460        ENDDO        ENDDO
# Line 283  C--   Upper diagonal Line 462  C--   Upper diagonal
462        DO j=jMin,jMax        DO j=jMin,jMax
463         DO i=iMin,iMax         DO i=iMin,iMax
464           c(i,j,1)  = 0. _d 0           c(i,j,1)  = 0. _d 0
          c(i,j,Nr) = 0. _d 0  
465         ENDDO         ENDDO
466        ENDDO        ENDDO
467  CML      DO k=1,Nr-1        DO k=2,Nr
       DO k=2,Nr-1  
        kp1=min(Nr,k+1)  
468         DO j=jMin,jMax         DO j=jMin,jMax
469          DO i=iMin,iMax          DO i=iMin,iMax
470              kp1=MAX(1,MIN(klowC(i,j,bi,bj),k+1))
471    C-    We keep recip_hFacC in the diffusive flux calculation,
472    C-    but no hFacC in TKE volume control
473    C-    No need for maskC(k) with recip_hFacC(k)
474            c(i,j,k) = -deltaTggl90            c(i,j,k) = -deltaTggl90
475       &        *recip_drF( k )*recip_hFacI(i,j,k,bi,bj)       &        *recip_drF( k ) * recip_hFacC(i,j,k,bi,bj)
476       &               *.5*(KappaE(i,j,k)+KappaE(i,j,kp1))       &        *.5 _d 0*(KappaE(i,j,k)+KappaE(i,j,kp1))
477       &        *recip_drC(k)       &        *recip_drC(k)*maskC(i,j,k-1,bi,bj)
           IF (recip_hFacI(i,j,kp1,bi,bj).EQ.0.) c(i,j,k)=0.  
478          ENDDO          ENDDO
479         ENDDO         ENDDO
480        ENDDO        ENDDO
481  C--   Center diagonal  C--   Center diagonal
482        DO k=1,Nr        DO k=1,Nr
483           km1 = MAX(k-1,1)
484         DO j=jMin,jMax         DO j=jMin,jMax
485          DO i=iMin,iMax          DO i=iMin,iMax
486            b(i,j,k) = 1. _d 0 - c(i,j,k) - a(i,j,k)            b(i,j,k) = 1. _d 0 - c(i,j,k) - a(i,j,k)
487       &        + ab15*deltaTggl90*GGL90ceps*SQRT(GGL90TKE(I,J,K,bi,bj))       &        + ab15*deltaTggl90*GGL90ceps*SQRTTKE(I,J,K)
488       &        /GGL90mixingLength(I,J,K)       &        * rMixingLength(I,J,K)
489         &        * maskC(i,j,k,bi,bj)*maskC(i,j,km1,bi,bj)
490           ENDDO           ENDDO
491         ENDDO         ENDDO
492        ENDDO        ENDDO
493  C     end set up matrix  C     end set up matrix
494    
 C  
495  C     Apply boundary condition  C     Apply boundary condition
496  C            kp1 = MIN(Nr,kSurf+1)
497        DO J=jMin,jMax        DO J=jMin,jMax
498         DO I=iMin,iMax         DO I=iMin,iMax
499  C     estimate friction velocity uStar from surface forcing  C     estimate friction velocity uStar from surface forcing
500          uStarSquare = SQRT(          uStarSquare = SQRT(
501       &         ( .5*( surfaceForcingU(I,  J,  bi,bj)       &    ( .5 _d 0*( surfaceForcingU(I,  J,  bi,bj)
502       &              + surfaceForcingU(I+1,J,  bi,bj) ) )**2       &              + surfaceForcingU(I+1,J,  bi,bj) ) )**2
503       &       + ( .5*( surfaceForcingV(I,  J,  bi,bj)       &  + ( .5 _d 0*( surfaceForcingV(I,  J,  bi,bj)
504       &              + surfaceForcingV(I,  J+1,bi,bj) ) )**2       &              + surfaceForcingV(I,  J+1,bi,bj) ) )**2
505       &                     )       &                     )
506  C     Dirichlet surface boundary condition for TKE  C     Dirichlet surface boundary condition for TKE
507          gTKE(I,J,kSurf) = MAX(GGL90TKEmin,GGL90m2*uStarSquare)          gTKE(I,J,kSurf) = MAX(GGL90TKEsurfMin,GGL90m2*uStarSquare)
508       &                     *maskC(I,J,kSurf,bi,bj)       &                     *maskC(I,J,kSurf,bi,bj)
509            gTKE(i,j,kp1) = gTKE(i,j,kp1)
510         &                - a(i,j,kp1)*gTKE(i,j,kSurf)
511            a(i,j,kp1) = 0. _d 0
512  C     Dirichlet bottom boundary condition for TKE = GGL90TKEbottom  C     Dirichlet bottom boundary condition for TKE = GGL90TKEbottom
513          kBottom   = MIN(MAX(kLowC(I,J,bi,bj),1),Nr)          kBottom   = MAX(kLowC(I,J,bi,bj),1)
514          gTKE(I,J,kBottom) = gTKE(I,J,kBottom)          gTKE(I,J,kBottom) = gTKE(I,J,kBottom)
515       &       - GGL90TKEbottom*c(I,J,kBottom)       &                    - GGL90TKEbottom*c(I,J,kBottom)
516            c(I,J,kBottom) = 0. _d 0
517         ENDDO         ENDDO
518        ENDDO            ENDDO
519  C  
520  C     solve tri-diagonal system, and store solution on gTKE (previously rhs)  C     solve tri-diagonal system, and store solution on gTKE (previously rhs)
521  C        CALL SOLVE_TRIDIAGONAL( iMin,iMax, jMin,jMax,
522        CALL GGL90_SOLVE( bi, bj, iMin, iMax, jMin, jMax,       I                        a, b, c,
523       I     a, b, c,       U                        gTKE,
524       U     gTKE,       O                        errCode,
525       I     myThid )       I                        bi, bj, myThid )
526  C  
527  C     now update TKE  C     now update TKE
 C      
528        DO K=1,Nr        DO K=1,Nr
529         DO J=jMin,jMax         DO J=jMin,jMax
530          DO I=iMin,iMax          DO I=iMin,iMax
531  C     impose minimum TKE to avoid numerical undershoots below zero  C     impose minimum TKE to avoid numerical undershoots below zero
532           GGL90TKE(I,J,K,bi,bj) = MAX( gTKE(I,J,K), GGL90TKEmin )           GGL90TKE(I,J,K,bi,bj) = MAX( gTKE(I,J,K), GGL90TKEmin )
533       &        * maskC(I,J,K,bi,bj)       &        * maskC(I,J,K,bi,bj)
534          ENDDO          ENDDO
535         ENDDO         ENDDO
536        ENDDO            ENDDO
537  C  
538  C     end of time step  C     end of time step
539  C     ===============================  C     ===============================
540  C     compute viscosity coefficients  
541  C            DO K=1,Nr
       DO K=2,Nr  
542         DO J=jMin,jMax         DO J=jMin,jMax
543          DO I=iMin,iMax          DO I=iMin,iMax
544  C     Eq. (11), (18) and (21)  #ifdef ALLOW_GGL90_SMOOTH
545           KappaM = GGL90ck*GGL90mixingLength(I,J,K)*           tmpVisc=
546       &                  SQRT( GGL90TKE(I,J,K,bi,bj) )       &  (
547           KappaH = KappaM/TKEPrandtlNumber(I,J,K)       &   p4 *  GGL90visctmp(i  ,j  ,k) * mskCor(i  ,j  ,bi,bj)
548  C     Set a minium (= background) value       &  +p8 *( GGL90visctmp(i-1,j  ,k) * mskCor(i-1,j  ,bi,bj)
549           KappaM = MAX(KappaM,viscAr)       &       + GGL90visctmp(i  ,j-1,k) * mskCor(i  ,j-1,bi,bj)
550           KappaH = MAX(KappaH,diffKrNrT(k))       &       + GGL90visctmp(i+1,j  ,k) * mskCor(i+1,j  ,bi,bj)
551  C     Set a maximum and mask land point       &       + GGL90visctmp(i  ,j+1,k) * mskCor(i  ,j+1,bi,bj))
552           GGL90viscAr(I,J,K,bi,bj) = MIN(KappaM,GGL90viscMax)       &  +p16*( GGL90visctmp(i+1,j+1,k) * mskCor(i+1,j+1,bi,bj)
553       &        * maskC(I,J,K,bi,bj)       &       + GGL90visctmp(i+1,j-1,k) * mskCor(i+1,j-1,bi,bj)
554           GGL90diffKr(I,J,K,bi,bj) = MIN(KappaH,GGL90diffMax)       &       + GGL90visctmp(i-1,j+1,k) * mskCor(i-1,j+1,bi,bj)
555       &        * maskC(I,J,K,bi,bj)       &       + GGL90visctmp(i-1,j-1,k) * mskCor(i-1,j-1,bi,bj))
556         &  )
557         & /(p4
558         &  +p8 *(       maskC(i-1,j  ,k,bi,bj) * mskCor(i-1,j  ,bi,bj)
559         &       +       maskC(i  ,j-1,k,bi,bj) * mskCor(i  ,j-1,bi,bj)
560         &       +       maskC(i+1,j  ,k,bi,bj) * mskCor(i+1,j  ,bi,bj)
561         &       +       maskC(i  ,j+1,k,bi,bj) * mskCor(i  ,j+1,bi,bj))
562         &  +p16*(       maskC(i+1,j+1,k,bi,bj) * mskCor(i+1,j+1,bi,bj)
563         &       +       maskC(i+1,j-1,k,bi,bj) * mskCor(i+1,j-1,bi,bj)
564         &       +       maskC(i-1,j+1,k,bi,bj) * mskCor(i-1,j+1,bi,bj)
565         &       +       maskC(i-1,j-1,k,bi,bj) * mskCor(i-1,j-1,bi,bj))
566         &  )*maskC(i,j,k,bi,bj)*mskCor(i,j,bi,bj)
567    #else
568             tmpVisc = GGL90visctmp(I,J,K)
569    #endif
570             tmpVisc = MIN(tmpVisc/TKEPrandtlNumber(i,j,k),GGL90diffMax)
571             GGL90diffKr(I,J,K,bi,bj)= MAX( tmpVisc , diffKrNrT(k) )
572          ENDDO          ENDDO
573         ENDDO         ENDDO
574  C     end third k-loop        ENDDO
575        ENDDO      
576    
577    
578          DO K=1,Nr
579           DO J=jMin,jMax
580            DO I=iMin,iMax
581    #ifdef ALLOW_GGL90_SMOOTH
582            tmpVisc =
583         & (
584         &   p4 *(GGL90visctmp(i  ,j  ,k) * mskCor(i  ,j  ,bi,bj)
585         &       +GGL90visctmp(i-1,j  ,k) * mskCor(i-1,j  ,bi,bj))
586         &  +p8 *(GGL90visctmp(i-1,j-1,k) * mskCor(i-1,j-1,bi,bj)
587         &       +GGL90visctmp(i-1,j+1,k) * mskCor(i-1,j+1,bi,bj)
588         &       +GGL90visctmp(i  ,j-1,k) * mskCor(i  ,j-1,bi,bj)
589         &       +GGL90visctmp(i  ,j+1,k) * mskCor(i  ,j+1,bi,bj))
590         &  )
591         & /(p4 * 2. _d 0
592         &  +p8 *(      maskC(i-1,j-1,k,bi,bj) * mskCor(i-1,j-1,bi,bj)
593         &       +      maskC(i-1,j+1,k,bi,bj) * mskCor(i-1,j+1,bi,bj)
594         &       +      maskC(i  ,j-1,k,bi,bj) * mskCor(i  ,j-1,bi,bj)
595         &       +      maskC(i  ,j+1,k,bi,bj) * mskCor(i  ,j+1,bi,bj))
596         &  )
597         &  *maskC(i  ,j,k,bi,bj)*mskCor(i  ,j,bi,bj)
598         &  *maskC(i-1,j,k,bi,bj)*mskCor(i-1,j,bi,bj)
599    #else
600            tmpVisc = _maskW(i,j,k,bi,bj) *
601         &                   (.5 _d 0*(GGL90visctmp(i,j,k)
602         &                            +GGL90visctmp(i-1,j,k))
603         &                   )
604    #endif
605            tmpVisc = MIN( tmpVisc , GGL90viscMax )
606            GGL90viscArU(i,j,k,bi,bj) = MAX( tmpVisc , viscArNr(k)  )
607            ENDDO
608           ENDDO
609          ENDDO
610    
611    
612          DO K=1,Nr
613           DO J=jMin,jMax
614            DO I=iMin,iMax
615    #ifdef ALLOW_GGL90_SMOOTH
616            tmpVisc =
617         & (
618         &   p4 *(GGL90visctmp(i  ,j  ,k) * mskCor(i  ,j  ,bi,bj)
619         &       +GGL90visctmp(i  ,j-1,k) * mskCor(i  ,j-1,bi,bj))
620         &  +p8 *(GGL90visctmp(i-1,j  ,k) * mskCor(i-1,j  ,bi,bj)
621         &       +GGL90visctmp(i-1,j-1,k) * mskCor(i-1,j-1,bi,bj)
622         &       +GGL90visctmp(i+1,j  ,k) * mskCor(i+1,j  ,bi,bj)
623         &       +GGL90visctmp(i+1,j-1,k) * mskCor(i+1,j-1,bi,bj))
624         &  )
625         & /(p4 * 2. _d 0
626         &  +p8 *(      maskC(i-1,j  ,k,bi,bj) * mskCor(i-1,j  ,bi,bj)
627         &       +      maskC(i-1,j-1,k,bi,bj) * mskCor(i-1,j-1,bi,bj)
628         &       +      maskC(i+1,j  ,k,bi,bj) * mskCor(i+1,j  ,bi,bj)
629         &       +      maskC(i+1,j-1,k,bi,bj) * mskCor(i+1,j-1,bi,bj))
630         &  )
631         &   *maskC(i,j  ,k,bi,bj)*mskCor(i,j  ,bi,bj)
632         &   *maskC(i,j-1,k,bi,bj)*mskCor(i,j-1,bi,bj)
633    #else
634            tmpVisc = _maskS(i,j,k,bi,bj) *
635         &                   (.5 _d 0*(GGL90visctmp(i,j,k)
636         &                            +GGL90visctmp(i,j-1,k))
637         &                   )
638    
639    #endif
640            tmpVisc = MIN( tmpVisc , GGL90viscMax )
641            GGL90viscArV(i,j,k,bi,bj) = MAX( tmpVisc , viscArNr(k)  )
642            ENDDO
643           ENDDO
644          ENDDO
645    
646    #ifdef ALLOW_DIAGNOSTICS
647          IF ( useDiagnostics ) THEN
648             CALL DIAGNOSTICS_FILL( GGL90TKE   ,'GGL90TKE',
649         &                          0,Nr, 1, bi, bj, myThid )
650             CALL DIAGNOSTICS_FILL( GGL90viscArU,'GGL90ArU',
651         &                          0,Nr, 1, bi, bj, myThid )
652             CALL DIAGNOSTICS_FILL( GGL90viscArV,'GGL90ArV',
653         &                          0,Nr, 1, bi, bj, myThid )
654             CALL DIAGNOSTICS_FILL( GGL90diffKr,'GGL90Kr ',
655         &                          0,Nr, 1, bi, bj, myThid )
656             CALL DIAGNOSTICS_FILL( TKEPrandtlNumber ,'GGL90Prl',
657         &                          0,Nr, 2, bi, bj, myThid )
658             CALL DIAGNOSTICS_FILL( GGL90mixingLength,'GGL90Lmx',
659         &                          0,Nr, 2, bi, bj, myThid )
660          ENDIF
661    #endif
662    
663  #endif /* ALLOW_GGL90 */  #endif /* ALLOW_GGL90 */
664    
665        RETURN        RETURN
666        END        END
   

Legend:
Removed from v.1.3  
changed lines
  Added in v.1.16

  ViewVC Help
Powered by ViewVC 1.1.22