/[MITgcm]/MITgcm/model/src/calc_oce_mxlayer.F
ViewVC logotype

Contents of /MITgcm/model/src/calc_oce_mxlayer.F

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


Revision 1.9 - (show annotations) (download)
Wed Jul 13 23:02:47 2011 UTC (12 years, 10 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint64i, checkpoint64h, checkpoint64a, checkpoint64c, checkpoint64b, checkpoint64e, checkpoint64d, checkpoint64g, checkpoint64f, checkpoint63p, checkpoint63q, checkpoint63r, checkpoint63s, checkpoint63l, checkpoint63m, checkpoint63n, checkpoint63o, checkpoint63h, checkpoint63i, checkpoint63j, checkpoint63k, checkpoint63d, checkpoint63e, checkpoint63f, checkpoint63g, checkpoint63a, checkpoint63b, checkpoint63c, checkpoint64
Changes since 1.8: +3 -3 lines
add Sub-Meso Eddies parameterisation (from Baylor)

1 C $Header: /u/gcmpack/MITgcm/model/src/calc_oce_mxlayer.F,v 1.8 2009/10/16 17:27:33 dimitri Exp $
2 C $Name: $
3
4 #include "PACKAGES_CONFIG.h"
5 #include "CPP_OPTIONS.h"
6
7 CBOP
8 C !ROUTINE: CALC_OCE_MXLAYER
9 C !INTERFACE:
10 SUBROUTINE CALC_OCE_MXLAYER(
11 I rhoSurf, sigmaR,
12 I bi, bj, myTime, myIter, myThid )
13 C !DESCRIPTION: \bv
14 C *==========================================================*
15 C | S/R CALC_OCE_MXLAYER
16 C | o Diagnose the Oceanic surface Mixed-Layer
17 C *==========================================================*
18 C \ev
19
20 C !USES:
21 IMPLICIT NONE
22 C == Global variables ==
23 #include "SIZE.h"
24 #include "EEPARAMS.h"
25 #include "PARAMS.h"
26 #include "GRID.h"
27 #include "DYNVARS.h"
28 #ifdef ALLOW_GMREDI
29 # include "GMREDI.h"
30 #endif
31
32 C !INPUT/OUTPUT PARAMETERS:
33 C == Routine Arguments ==
34 C rhoSurf :: Surface density anomaly
35 C sigmaR :: Vertical gradient of potential density
36 C bi,bj :: tile indices
37 C myTime :: Current time in simulation
38 C myIter :: Current iteration number in simulation
39 C myThid :: my Thread Id number
40 _RL rhoSurf(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
41 _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
42 INTEGER bi, bj
43 _RL myTime
44 INTEGER myIter
45 INTEGER myThid
46
47 C === Functions ====
48 #ifdef ALLOW_DIAGNOSTICS
49 LOGICAL DIAGNOSTICS_IS_ON
50 EXTERNAL DIAGNOSTICS_IS_ON
51 #endif /* ALLOW_DIAGNOSTICS */
52
53 C !LOCAL VARIABLES:
54 C == Local variables ==
55 C i,j :: Loop counters
56 INTEGER i,j,k
57 LOGICAL calcMixLayerDepth
58 INTEGER method
59 _RL rhoBigNb
60 _RL rhoMxL(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
61 _RL rhoKm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
62 _RL rhoLoc(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
63 _RL tmpFac, sigmAv
64 CEOP
65
66 calcMixLayerDepth = .FALSE.
67 #ifdef ALLOW_GMREDI
68 IF ( useGMRedi .AND. .NOT.useKPP ) THEN
69 calcMixLayerDepth = GM_useSubMeso .OR. GM_taper_scheme.EQ.'fm07'
70 ENDIF
71 #endif
72 #ifdef ALLOW_DIAGNOSTICS
73 IF ( useDiagnostics.AND. .NOT.calcMixLayerDepth ) THEN
74 calcMixLayerDepth = DIAGNOSTICS_IS_ON('MXLDEPTH',myThid)
75 ENDIF
76 #endif
77 IF ( calcMixLayerDepth ) THEN
78
79 C-- Select which "method" to use:
80 method = 0
81 IF ( hMixCriteria.LT.0. ) method = 1
82 IF ( hMixCriteria.GT.1. ) method = 2
83
84 IF ( method.EQ.1 ) THEN
85
86 C-- First method :
87 C where the potential density (ref.lev=surface) is larger than
88 C surface density plus Delta_Rho = hMixCriteria * Alpha(surf)
89 C = density of water which is -hMixCriteria colder than surface water
90 C (see Kara, Rochford, and Hurlburt JGR 2000 for default criterion)
91
92 c hMixCriteria = -0.8 _d 0
93 c dRhoSmall = 1. _d -6
94 rhoBigNb = rhoConst*1. _d 10
95 CALL FIND_ALPHA(
96 I bi, bj, 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1, 1,
97 O rhoMxL, myThid )
98
99 DO j=1-Oly,sNy+Oly
100 DO i=1-Olx,sNx+Olx
101 rhoKm1(i,j) = rhoSurf(i,j)
102 rhoMxL(i,j) = rhoSurf(i,j)
103 & + MAX( rhoMxL(i,j)*hMixCriteria, dRhoSmall )
104 hMixLayer(i,j,bi,bj) = rF(1)-R_low(i,j,bi,bj)
105 ENDDO
106 ENDDO
107 DO k = 2,Nr
108 C- potential density (reference level = surface level)
109 CALL FIND_RHO_2D(
110 I 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1,
111 I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
112 O rhoLoc,
113 I k, bi, bj, myThid )
114
115 DO j=1-Oly,sNy+Oly
116 DO i=1-Olx,sNx+Olx
117 IF ( k.LE.klowC(i,j,bi,bj) .AND.
118 & rhoLoc(i,j).GE.rhoMxL(i,j) ) THEN
119 IF ( rhoLoc(i,j).GT.rhoKm1(i,j) ) THEN
120 tmpFac = ( rhoMxL(i,j) - rhoKm1(i,j) )
121 & / ( rhoLoc(i,j) - rhoKm1(i,j) )
122 ELSE
123 tmpFac = 0.
124 ENDIF
125 hMixLayer(i,j,bi,bj) = rF(1)-rC(k-1)+tmpFac*drC(k)
126 rhoMxL(i,j) = rhoBigNb
127 ELSE
128 rhoKm1(i,j) = rhoLoc(i,j)
129 ENDIF
130 ENDDO
131 ENDDO
132 ENDDO
133
134 ELSEIF ( method.EQ.2 ) THEN
135
136 C-- Second method :
137 C where the local stratification exceed the mean stratification above
138 C (from surface down to here) by factor hMixCriteria
139
140 c hMixCriteria = 1.5 _d 0
141 c dRhoSmall = 1. _d -2
142 DO j=1-Oly,sNy+Oly
143 DO i=1-Olx,sNx+Olx
144 IF ( klowC(i,j,bi,bj) .GT. 0 ) THEN
145 hMixLayer(i,j,bi,bj) = drF(1)
146 rhoMxL(i,j) = 1.
147 ELSE
148 hMixLayer(i,j,bi,bj) = rF(1)
149 rhoMxL(i,j) = -1.
150 ENDIF
151 ENDDO
152 ENDDO
153 DO k = 2,Nr-1
154 C- potential density (reference level = surface level)
155 CALL FIND_RHO_2D(
156 I 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1,
157 I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
158 O rhoLoc,
159 I k, bi, bj, myThid )
160
161 DO j=1-Oly,sNy+Oly
162 DO i=1-Olx,sNx+Olx
163 IF ( k.LT.klowC(i,j,bi,bj)
164 & .AND. rhoMxL(i,j).GE.0. ) THEN
165 sigmAv = ( rhoLoc(i,j)-rhoSurf(i,j)+dRhoSmall )
166 & / ( rC(1)-rC(k) )
167 IF ( -sigmaR(i,j,k+1).GT.sigmAv*hMixCriteria ) THEN
168 tmpFac = 0. _d 0
169 IF ( sigmAv.GT.0. _d 0 ) THEN
170 tmpFac = hMixCriteria*sigmaR(i,j,k)/sigmaR(i,j,k+1)
171 IF ( tmpFac .GT. 1. _d 0 ) THEN
172 tmpFac = 1. _d 0
173 & + ( tmpFac - 1. _d 0 )/( hMixCriteria - 1. _d 0 )
174 ENDIF
175 tmpFac = MAX( 0. _d 0, MIN( tmpFac, 2. _d 0 ) )
176 ENDIF
177 hMixLayer(i,j,bi,bj) = rF(1)-rF(k+1)
178 & - drF(k)*tmpFac*0.5 _d 0
179 rhoMxL(i,j) = -1.
180 ENDIF
181 ENDIF
182 ENDDO
183 ENDDO
184 ENDDO
185
186 ELSE
187 STOP 'S/R CALC_OCE_MXLAYER: invalid method'
188 ENDIF
189
190 IF ( hMixSmooth .GT. 0. _d 0 ) THEN
191 tmpFac = (1. _d 0 - hMixSmooth ) / 4. _d 0
192 DO j=1-Oly+1,sNy+Oly-1
193 DO i=1-Olx+1,sNx+Olx-1
194 rhoLoc(i,j)=(hMixSmooth * hMixLayer(i,j,bi,bj) +
195 & tmpFac * ( hMixLayer(i-1,j,bi,bj) +
196 & hMixLayer(i+1,j,bi,bj) +
197 & hMixLayer(i,j-1,bi,bj) +
198 & hMixLayer(i,j+1,bi,bj) )
199 & )
200 & /(hMixSmooth +
201 & tmpFac * ( maskC(i-1,j,1,bi,bj) +
202 & maskC(i+1,j,1,bi,bj) +
203 & maskC(i,j-1,1,bi,bj) +
204 & maskC(i,j+1,1,bi,bj) )
205 & ) * maskC(i,j,1,bi,bj)
206 ENDDO
207 ENDDO
208 DO j=1-Oly+1,sNy+Oly-1
209 DO i=1-Olx+1,sNx+Olx-1
210 hMixLayer(i,j,bi,bj) = rhoLoc(i,j)
211 ENDDO
212 ENDDO
213 ENDIF
214
215 #ifdef ALLOW_DIAGNOSTICS
216 IF ( useDiagnostics ) THEN
217 CALL DIAGNOSTICS_FILL( hMixLayer, 'MXLDEPTH',
218 & 0, 1, 1, bi, bj, myThid )
219 ENDIF
220 #endif /* ALLOW_DIAGNOSTICS */
221
222 C-- end if calcMixLayerDepth
223 ENDIF
224
225 RETURN
226 END

  ViewVC Help
Powered by ViewVC 1.1.22