/[MITgcm]/MITgcm_contrib/submesoscale/code/calc_oce_mxlayer.F
ViewVC logotype

Diff of /MITgcm_contrib/submesoscale/code/calc_oce_mxlayer.F

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

revision 1.2 by dimitri, Fri May 30 23:18:08 2008 UTC revision 1.3 by zhc, Fri Mar 12 18:31:00 2010 UTC
# Line 24  C     == Global variables == Line 24  C     == Global variables ==
24  #include "EEPARAMS.h"  #include "EEPARAMS.h"
25  #include "PARAMS.h"  #include "PARAMS.h"
26  #include "GRID.h"  #include "GRID.h"
 #include "SURFACE.h"  
27  #include "DYNVARS.h"  #include "DYNVARS.h"
28  #ifdef ALLOW_GMREDI  #ifdef ALLOW_GMREDI
29  # include "GMREDI.h"  # include "GMREDI.h"
# Line 57  C     i,j :: Loop counters Line 56  C     i,j :: Loop counters
56        INTEGER i,j,k        INTEGER i,j,k
57        LOGICAL calcMixLayerDepth        LOGICAL calcMixLayerDepth
58        INTEGER method        INTEGER method
59        _RL     dRhoSmall, rhoBigNb        _RL     rhoBigNb
60        _RL     rhoMxL(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL     rhoMxL(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
61        _RL     rhoKm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL     rhoKm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
62        _RL     rhoLoc(1-OLx:sNx+OLx,1-OLy:sNy+OLy)        _RL     rhoLoc(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
# Line 66  CEOP Line 65  CEOP
65    
66        calcMixLayerDepth = .FALSE.        calcMixLayerDepth = .FALSE.
67  #ifdef ALLOW_GMREDI  #ifdef ALLOW_GMREDI
68        IF ( useGMRedi ) calcMixLayerDepth = GM_taper_scheme.EQ.'fm07'        IF ( useGMRedi ) THEN
69        IF ( useGMRedi .AND. GM_SM_Ce.GT.0.0 ) calcMixLayerDepth = .TRUE.         IF ( .NOT.useKPP ) calcMixLayerDepth = GM_taper_scheme.EQ.'fm07'
70    C SM(1)
71          IF ( GM_SM_Ce.GT.0.0 ) calcMixLayerDepth = .TRUE.
72          ENDIF
73  #endif  #endif
74  #ifdef ALLOW_DIAGNOSTICS  #ifdef ALLOW_DIAGNOSTICS
75        IF ( useDiagnostics.AND. .NOT.calcMixLayerDepth ) THEN        IF ( useDiagnostics.AND. .NOT.calcMixLayerDepth ) THEN
# Line 87  C--   First method : Line 89  C--   First method :
89  C     where the potential density (ref.lev=surface) is larger than  C     where the potential density (ref.lev=surface) is larger than
90  C       surface density plus Delta_Rho = hMixCriteria * Alpha(surf)  C       surface density plus Delta_Rho = hMixCriteria * Alpha(surf)
91  C     = density of water which is -hMixCriteria colder than surface water  C     = density of water which is -hMixCriteria colder than surface water
92    C     (see Kara, Rochford, and Hurlburt JGR 2000 for default criterion)
93    
94  c       hMixCriteria  = -0.8 _d 0  c       hMixCriteria  = -0.8 _d 0
95          dRhoSmall = 1. _d -6  c       dRhoSmall = 1. _d -6
96          rhoBigNb  = rhoConst*1. _d 10          rhoBigNb  = rhoConst*1. _d 10
97          CALL FIND_ALPHA(          CALL FIND_ALPHA(
98       I            bi, bj, 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1, 1,       I            bi, bj, 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1, 1,
# Line 100  c       hMixCriteria  = -0.8 _d 0 Line 103  c       hMixCriteria  = -0.8 _d 0
103             rhoKm1(i,j) = rhoSurf(i,j)             rhoKm1(i,j) = rhoSurf(i,j)
104             rhoMxL(i,j) = rhoSurf(i,j)             rhoMxL(i,j) = rhoSurf(i,j)
105       &                 + MAX( rhoMxL(i,j)*hMixCriteria, dRhoSmall )       &                 + MAX( rhoMxL(i,j)*hMixCriteria, dRhoSmall )
106             hMixLayer(i,j,bi,bj) = rF(1)-R_low(I,J,bi,bj)             hMixLayer(i,j,bi,bj) = rF(1)-R_low(i,j,bi,bj)
107           ENDDO           ENDDO
108          ENDDO          ENDDO
109          DO k = 2,Nr          DO k = 2,Nr
110  C-    potential density (reference level = surface level)  C-    potential density (reference level = surface level)
111           CALL FIND_RHO(           CALL FIND_RHO_2D(
112       I             bi, bj, 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, K, 1,       I        1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1,
113       I             theta, salt,       I        theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
114       O             rhoLoc, myThid )       O        rhoLoc,
115         I        k, bi, bj, myThid )
116    
117           DO j=1-Oly,sNy+Oly           DO j=1-Oly,sNy+Oly
118            DO i=1-Olx,sNx+Olx            DO i=1-Olx,sNx+Olx
# Line 135  C--   Second method : Line 139  C--   Second method :
139  C     where the local stratification exceed the mean stratification above  C     where the local stratification exceed the mean stratification above
140  C     (from surface down to here) by factor hMixCriteria  C     (from surface down to here) by factor hMixCriteria
141    
142  c       hMixCriteria  = 2. _d 0  c       hMixCriteria  = 1.5 _d 0
143  C-Note: dRhoSmall is hard coded for now but should become run-time parameter  c       dRhoSmall = 1. _d -2
         dRhoSmall = 1. _d -2  
144          DO j=1-Oly,sNy+Oly          DO j=1-Oly,sNy+Oly
145           DO i=1-Olx,sNx+Olx           DO i=1-Olx,sNx+Olx
146             IF ( klowC(i,j,bi,bj) .GT. 0 ) THEN             IF ( klowC(i,j,bi,bj) .GT. 0 ) THEN
# Line 151  C-Note: dRhoSmall is hard coded for now Line 154  C-Note: dRhoSmall is hard coded for now
154          ENDDO          ENDDO
155          DO k = 2,Nr-1          DO k = 2,Nr-1
156  C-    potential density (reference level = surface level)  C-    potential density (reference level = surface level)
157           CALL FIND_RHO(           CALL FIND_RHO_2D(
158       I             bi, bj, 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, K, 1,       I        1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1,
159       I             theta, salt,       I        theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
160       O             rhoLoc, myThid )       O        rhoLoc,
161         I        k, bi, bj, myThid )
162    
163           DO j=1-Oly,sNy+Oly           DO j=1-Oly,sNy+Oly
164            DO i=1-Olx,sNx+Olx            DO i=1-Olx,sNx+Olx
# Line 185  C-    potential density (reference level Line 189  C-    potential density (reference level
189          STOP 'S/R CALC_OCE_MXLAYER: invalid method'          STOP 'S/R CALC_OCE_MXLAYER: invalid method'
190         ENDIF         ENDIF
191    
192           IF ( hMixSmooth .GT. 0. _d 0 ) THEN
193            tmpFac = (1. _d 0 - hMixSmooth ) / 4. _d 0
194            DO j=1-Oly+1,sNy+Oly-1
195             DO i=1-Olx+1,sNx+Olx-1
196                rhoLoc(i,j)=(hMixSmooth *   hMixLayer(i,j,bi,bj)   +
197         &                       tmpFac * ( hMixLayer(i-1,j,bi,bj) +
198         &                                  hMixLayer(i+1,j,bi,bj) +
199         &                                  hMixLayer(i,j-1,bi,bj) +
200         &                                  hMixLayer(i,j+1,bi,bj) )
201         &                  )
202         &                 /(hMixSmooth +
203         &                       tmpFac * ( maskC(i-1,j,1,bi,bj) +
204         &                                  maskC(i+1,j,1,bi,bj) +
205         &                                  maskC(i,j-1,1,bi,bj) +
206         &                                  maskC(i,j+1,1,bi,bj) )
207         &                  ) * maskC(i,j,1,bi,bj)
208             ENDDO
209            ENDDO
210            DO j=1-Oly+1,sNy+Oly-1
211             DO i=1-Olx+1,sNx+Olx-1
212                hMixLayer(i,j,bi,bj) = rhoLoc(i,j)
213             ENDDO
214            ENDDO
215           ENDIF
216    
217  #ifdef ALLOW_DIAGNOSTICS  #ifdef ALLOW_DIAGNOSTICS
218         IF ( useDiagnostics ) THEN         IF ( useDiagnostics ) THEN
219          CALL DIAGNOSTICS_FILL( hMixLayer, 'MXLDEPTH',          CALL DIAGNOSTICS_FILL( hMixLayer, 'MXLDEPTH',

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

  ViewVC Help
Powered by ViewVC 1.1.22