/[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.3 - (show annotations) (download)
Mon Aug 11 22:25:52 2008 UTC (15 years, 9 months ago) by jmc
Branch: MAIN
Changes since 1.2: +13 -11 lines
replace calls to "FIND_RHO" with recent version "FIND_RHO_2D"

1 C $Header: /u/gcmpack/MITgcm/model/src/calc_oce_mxlayer.F,v 1.2 2007/06/20 23:06:51 jmc 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 "SURFACE.h"
28 #include "DYNVARS.h"
29 #ifdef ALLOW_GMREDI
30 # include "GMREDI.h"
31 #endif
32
33 C !INPUT/OUTPUT PARAMETERS:
34 C == Routine Arguments ==
35 C rhoSurf :: Surface density anomaly
36 C sigmaR :: Vertical gradient of potential density
37 C bi,bj :: tile indices
38 C myTime :: Current time in simulation
39 C myIter :: Current iteration number in simulation
40 C myThid :: my Thread Id number
41 _RL rhoSurf(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
42 _RL sigmaR (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
43 INTEGER bi, bj
44 _RL myTime
45 INTEGER myIter
46 INTEGER myThid
47
48 C === Functions ====
49 #ifdef ALLOW_DIAGNOSTICS
50 LOGICAL DIAGNOSTICS_IS_ON
51 EXTERNAL DIAGNOSTICS_IS_ON
52 #endif /* ALLOW_DIAGNOSTICS */
53
54 C !LOCAL VARIABLES:
55 C == Local variables ==
56 C i,j :: Loop counters
57 INTEGER i,j,k
58 LOGICAL calcMixLayerDepth
59 INTEGER method
60 _RL dRhoSmall, rhoBigNb
61 _RL rhoMxL(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
62 _RL rhoKm1(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
63 _RL rhoLoc(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
64 _RL tmpFac, sigmAv
65 CEOP
66
67 calcMixLayerDepth = .FALSE.
68 #ifdef ALLOW_GMREDI
69 IF ( useGMRedi ) calcMixLayerDepth = GM_taper_scheme.EQ.'fm07'
70 #endif
71 #ifdef ALLOW_DIAGNOSTICS
72 IF ( useDiagnostics.AND. .NOT.calcMixLayerDepth ) THEN
73 calcMixLayerDepth = DIAGNOSTICS_IS_ON('MXLDEPTH',myThid)
74 ENDIF
75 #endif
76 IF ( calcMixLayerDepth ) THEN
77
78 C-- Select which "method" to use:
79 method = 0
80 IF ( hMixCriteria.LT.0. ) method = 1
81 IF ( hMixCriteria.GT.1. ) method = 2
82
83 IF ( method.EQ.1 ) THEN
84
85 C-- First method :
86 C where the potential density (ref.lev=surface) is larger than
87 C surface density plus Delta_Rho = hMixCriteria * Alpha(surf)
88 C = density of water which is -hMixCriteria colder than surface water
89
90 c hMixCriteria = -0.8 _d 0
91 dRhoSmall = 1. _d -6
92 rhoBigNb = rhoConst*1. _d 10
93 CALL FIND_ALPHA(
94 I bi, bj, 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1, 1,
95 O rhoMxL, myThid )
96
97 DO j=1-Oly,sNy+Oly
98 DO i=1-Olx,sNx+Olx
99 rhoKm1(i,j) = rhoSurf(i,j)
100 rhoMxL(i,j) = rhoSurf(i,j)
101 & + MAX( rhoMxL(i,j)*hMixCriteria, dRhoSmall )
102 hMixLayer(i,j,bi,bj) = rF(1)-R_low(i,j,bi,bj)
103 ENDDO
104 ENDDO
105 DO k = 2,Nr
106 C- potential density (reference level = surface level)
107 CALL FIND_RHO_2D(
108 I 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1,
109 I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
110 O rhoLoc,
111 I k, bi, bj, myThid )
112
113 DO j=1-Oly,sNy+Oly
114 DO i=1-Olx,sNx+Olx
115 IF ( k.LE.klowC(i,j,bi,bj) .AND.
116 & rhoLoc(i,j).GE.rhoMxL(i,j) ) THEN
117 IF ( rhoLoc(i,j).GT.rhoKm1(i,j) ) THEN
118 tmpFac = ( rhoMxL(i,j) - rhoKm1(i,j) )
119 & / ( rhoLoc(i,j) - rhoKm1(i,j) )
120 ELSE
121 tmpFac = 0.
122 ENDIF
123 hMixLayer(i,j,bi,bj) = rF(1)-rC(k-1)+tmpFac*drC(k)
124 rhoMxL(i,j) = rhoBigNb
125 ELSE
126 rhoKm1(i,j) = rhoLoc(i,j)
127 ENDIF
128 ENDDO
129 ENDDO
130 ENDDO
131
132 ELSEIF ( method.EQ.2 ) THEN
133
134 C-- Second method :
135 C where the local stratification exceed the mean stratification above
136 C (from surface down to here) by factor hMixCriteria
137
138 c hMixCriteria = 2. _d 0
139 C-Note: dRhoSmall is hard coded for now but should become run-time parameter
140 dRhoSmall = 1. _d -2
141 DO j=1-Oly,sNy+Oly
142 DO i=1-Olx,sNx+Olx
143 IF ( klowC(i,j,bi,bj) .GT. 0 ) THEN
144 hMixLayer(i,j,bi,bj) = drF(1)
145 rhoMxL(i,j) = 1.
146 ELSE
147 hMixLayer(i,j,bi,bj) = rF(1)
148 rhoMxL(i,j) = -1.
149 ENDIF
150 ENDDO
151 ENDDO
152 DO k = 2,Nr-1
153 C- potential density (reference level = surface level)
154 CALL FIND_RHO_2D(
155 I 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1,
156 I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
157 O rhoLoc,
158 I k, bi, bj, myThid )
159
160 DO j=1-Oly,sNy+Oly
161 DO i=1-Olx,sNx+Olx
162 IF ( k.LT.klowC(i,j,bi,bj)
163 & .AND. rhoMxL(i,j).GE.0. ) THEN
164 sigmAv = ( rhoLoc(i,j)-rhoSurf(i,j)+dRhoSmall )
165 & / ( rC(1)-rC(k) )
166 IF ( -sigmaR(i,j,k+1).GT.sigmAv*hMixCriteria ) THEN
167 tmpFac = 0. _d 0
168 IF ( sigmAv.GT.0. _d 0 ) THEN
169 tmpFac = hMixCriteria*sigmaR(i,j,k)/sigmaR(i,j,k+1)
170 IF ( tmpFac .GT. 1. _d 0 ) THEN
171 tmpFac = 1. _d 0
172 & + ( tmpFac - 1. _d 0 )/( hMixCriteria - 1. _d 0 )
173 ENDIF
174 tmpFac = MAX( 0. _d 0, MIN( tmpFac, 2. _d 0 ) )
175 ENDIF
176 hMixLayer(i,j,bi,bj) = rF(1)-rF(k+1)
177 & - drF(k)*tmpFac*0.5 _d 0
178 rhoMxL(i,j) = -1.
179 ENDIF
180 ENDIF
181 ENDDO
182 ENDDO
183 ENDDO
184
185 ELSE
186 STOP 'S/R CALC_OCE_MXLAYER: invalid method'
187 ENDIF
188
189 #ifdef ALLOW_DIAGNOSTICS
190 IF ( useDiagnostics ) THEN
191 CALL DIAGNOSTICS_FILL( hMixLayer, 'MXLDEPTH',
192 & 0, 1, 1, bi, bj, myThid )
193 ENDIF
194 #endif /* ALLOW_DIAGNOSTICS */
195
196 C-- end if calcMixLayerDepth
197 ENDIF
198
199 RETURN
200 END

  ViewVC Help
Powered by ViewVC 1.1.22