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

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

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


Revision 1.10 - (hide annotations) (download)
Thu Jun 27 14:47:01 2013 UTC (10 years, 10 months ago) by m_bates
Branch: MAIN
CVS Tags: checkpoint64y, checkpoint64x, checkpoint64z, checkpoint64q, checkpoint64p, checkpoint64s, checkpoint64r, checkpoint64u, checkpoint64t, checkpoint64w, checkpoint64v, checkpoint64k, checkpoint64j, checkpoint64m, checkpoint64l, checkpoint64o, checkpoint64n, checkpoint65, checkpoint66g, checkpoint66f, checkpoint66e, checkpoint66d, checkpoint66c, checkpoint66b, checkpoint66a, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, checkpoint65z, checkpoint65x, checkpoint65y, checkpoint65r, checkpoint65s, checkpoint65p, checkpoint65q, checkpoint65v, checkpoint65w, checkpoint65t, checkpoint65u, checkpoint65j, checkpoint65k, checkpoint65h, checkpoint65i, checkpoint65n, checkpoint65o, checkpoint65l, checkpoint65m, checkpoint65b, checkpoint65c, checkpoint65a, checkpoint65f, checkpoint65g, checkpoint65d, checkpoint65e, HEAD
Changes since 1.9: +2 -1 lines
use calc_oce_mxlayer to calculate the mixed layer depth for the PV eddy closure in gmredi_k3d

1 m_bates 1.10 C $Header: /u/gcmpack/MITgcm/model/src/calc_oce_mxlayer.F,v 1.9 2011/07/13 23:02:47 jmc Exp $
2 jmc 1.1 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 dfer 1.6 _RL rhoBigNb
60 jmc 1.1 _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 jmc 1.9 IF ( useGMRedi .AND. .NOT.useKPP ) THEN
69     calcMixLayerDepth = GM_useSubMeso .OR. GM_taper_scheme.EQ.'fm07'
70 m_bates 1.10 & .OR. GM_useK3D
71 dfer 1.5 ENDIF
72 jmc 1.1 #endif
73     #ifdef ALLOW_DIAGNOSTICS
74     IF ( useDiagnostics.AND. .NOT.calcMixLayerDepth ) THEN
75     calcMixLayerDepth = DIAGNOSTICS_IS_ON('MXLDEPTH',myThid)
76     ENDIF
77     #endif
78     IF ( calcMixLayerDepth ) THEN
79    
80     C-- Select which "method" to use:
81     method = 0
82     IF ( hMixCriteria.LT.0. ) method = 1
83     IF ( hMixCriteria.GT.1. ) method = 2
84    
85     IF ( method.EQ.1 ) THEN
86    
87     C-- First method :
88     C where the potential density (ref.lev=surface) is larger than
89     C surface density plus Delta_Rho = hMixCriteria * Alpha(surf)
90     C = density of water which is -hMixCriteria colder than surface water
91 dimitri 1.8 C (see Kara, Rochford, and Hurlburt JGR 2000 for default criterion)
92 jmc 1.1
93     c hMixCriteria = -0.8 _d 0
94 dfer 1.6 c dRhoSmall = 1. _d -6
95 jmc 1.1 rhoBigNb = rhoConst*1. _d 10
96     CALL FIND_ALPHA(
97     I bi, bj, 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1, 1,
98     O rhoMxL, myThid )
99    
100     DO j=1-Oly,sNy+Oly
101     DO i=1-Olx,sNx+Olx
102     rhoKm1(i,j) = rhoSurf(i,j)
103     rhoMxL(i,j) = rhoSurf(i,j)
104     & + MAX( rhoMxL(i,j)*hMixCriteria, dRhoSmall )
105 jmc 1.3 hMixLayer(i,j,bi,bj) = rF(1)-R_low(i,j,bi,bj)
106 jmc 1.1 ENDDO
107     ENDDO
108     DO k = 2,Nr
109     C- potential density (reference level = surface level)
110 jmc 1.3 CALL FIND_RHO_2D(
111     I 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1,
112     I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
113     O rhoLoc,
114     I k, bi, bj, myThid )
115 jmc 1.1
116     DO j=1-Oly,sNy+Oly
117     DO i=1-Olx,sNx+Olx
118     IF ( k.LE.klowC(i,j,bi,bj) .AND.
119     & rhoLoc(i,j).GE.rhoMxL(i,j) ) THEN
120     IF ( rhoLoc(i,j).GT.rhoKm1(i,j) ) THEN
121     tmpFac = ( rhoMxL(i,j) - rhoKm1(i,j) )
122     & / ( rhoLoc(i,j) - rhoKm1(i,j) )
123     ELSE
124     tmpFac = 0.
125     ENDIF
126     hMixLayer(i,j,bi,bj) = rF(1)-rC(k-1)+tmpFac*drC(k)
127     rhoMxL(i,j) = rhoBigNb
128     ELSE
129     rhoKm1(i,j) = rhoLoc(i,j)
130     ENDIF
131     ENDDO
132     ENDDO
133     ENDDO
134    
135     ELSEIF ( method.EQ.2 ) THEN
136    
137     C-- Second method :
138     C where the local stratification exceed the mean stratification above
139     C (from surface down to here) by factor hMixCriteria
140    
141 dfer 1.6 c hMixCriteria = 1.5 _d 0
142     c dRhoSmall = 1. _d -2
143 jmc 1.1 DO j=1-Oly,sNy+Oly
144     DO i=1-Olx,sNx+Olx
145     IF ( klowC(i,j,bi,bj) .GT. 0 ) THEN
146     hMixLayer(i,j,bi,bj) = drF(1)
147     rhoMxL(i,j) = 1.
148     ELSE
149     hMixLayer(i,j,bi,bj) = rF(1)
150     rhoMxL(i,j) = -1.
151     ENDIF
152     ENDDO
153     ENDDO
154 jmc 1.2 DO k = 2,Nr-1
155 jmc 1.1 C- potential density (reference level = surface level)
156 jmc 1.3 CALL FIND_RHO_2D(
157     I 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1,
158     I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
159     O rhoLoc,
160     I k, bi, bj, myThid )
161 jmc 1.1
162     DO j=1-Oly,sNy+Oly
163     DO i=1-Olx,sNx+Olx
164 jmc 1.2 IF ( k.LT.klowC(i,j,bi,bj)
165 jmc 1.1 & .AND. rhoMxL(i,j).GE.0. ) THEN
166 jmc 1.2 sigmAv = ( rhoLoc(i,j)-rhoSurf(i,j)+dRhoSmall )
167     & / ( rC(1)-rC(k) )
168     IF ( -sigmaR(i,j,k+1).GT.sigmAv*hMixCriteria ) THEN
169 jmc 1.1 tmpFac = 0. _d 0
170     IF ( sigmAv.GT.0. _d 0 ) THEN
171 jmc 1.2 tmpFac = hMixCriteria*sigmaR(i,j,k)/sigmaR(i,j,k+1)
172     IF ( tmpFac .GT. 1. _d 0 ) THEN
173     tmpFac = 1. _d 0
174     & + ( tmpFac - 1. _d 0 )/( hMixCriteria - 1. _d 0 )
175     ENDIF
176     tmpFac = MAX( 0. _d 0, MIN( tmpFac, 2. _d 0 ) )
177 jmc 1.1 ENDIF
178 jmc 1.2 hMixLayer(i,j,bi,bj) = rF(1)-rF(k+1)
179     & - drF(k)*tmpFac*0.5 _d 0
180 jmc 1.1 rhoMxL(i,j) = -1.
181     ENDIF
182     ENDIF
183     ENDDO
184     ENDDO
185     ENDDO
186    
187     ELSE
188     STOP 'S/R CALC_OCE_MXLAYER: invalid method'
189     ENDIF
190    
191 dfer 1.7 IF ( hMixSmooth .GT. 0. _d 0 ) THEN
192     tmpFac = (1. _d 0 - hMixSmooth ) / 4. _d 0
193     DO j=1-Oly+1,sNy+Oly-1
194     DO i=1-Olx+1,sNx+Olx-1
195     rhoLoc(i,j)=(hMixSmooth * hMixLayer(i,j,bi,bj) +
196     & tmpFac * ( hMixLayer(i-1,j,bi,bj) +
197     & hMixLayer(i+1,j,bi,bj) +
198     & hMixLayer(i,j-1,bi,bj) +
199     & hMixLayer(i,j+1,bi,bj) )
200     & )
201     & /(hMixSmooth +
202     & tmpFac * ( maskC(i-1,j,1,bi,bj) +
203     & maskC(i+1,j,1,bi,bj) +
204     & maskC(i,j-1,1,bi,bj) +
205     & maskC(i,j+1,1,bi,bj) )
206     & ) * maskC(i,j,1,bi,bj)
207     ENDDO
208     ENDDO
209     DO j=1-Oly+1,sNy+Oly-1
210     DO i=1-Olx+1,sNx+Olx-1
211     hMixLayer(i,j,bi,bj) = rhoLoc(i,j)
212     ENDDO
213     ENDDO
214     ENDIF
215    
216 jmc 1.1 #ifdef ALLOW_DIAGNOSTICS
217     IF ( useDiagnostics ) THEN
218     CALL DIAGNOSTICS_FILL( hMixLayer, 'MXLDEPTH',
219     & 0, 1, 1, bi, bj, myThid )
220     ENDIF
221     #endif /* ALLOW_DIAGNOSTICS */
222    
223     C-- end if calcMixLayerDepth
224     ENDIF
225    
226     RETURN
227     END

  ViewVC Help
Powered by ViewVC 1.1.22