/[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.5 - (hide annotations) (download)
Wed Oct 29 22:28:51 2008 UTC (15 years, 6 months ago) by dfer
Branch: MAIN
Changes since 1.4: +4 -2 lines
Avoid computation if KPP turned on

1 dfer 1.5 C $Header: /u/gcmpack/MITgcm/model/src/calc_oce_mxlayer.F,v 1.4 2008/08/12 22:26:39 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     _RL dRhoSmall, 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 dfer 1.5 IF ( useGMRedi ) THEN
69     IF ( .NOT.useKPP ) calcMixLayerDepth = GM_taper_scheme.EQ.'fm07'
70     ENDIF
71 jmc 1.1 #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    
91     c hMixCriteria = -0.8 _d 0
92     dRhoSmall = 1. _d -6
93     rhoBigNb = rhoConst*1. _d 10
94     CALL FIND_ALPHA(
95     I bi, bj, 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1, 1,
96     O rhoMxL, myThid )
97    
98     DO j=1-Oly,sNy+Oly
99     DO i=1-Olx,sNx+Olx
100     rhoKm1(i,j) = rhoSurf(i,j)
101     rhoMxL(i,j) = rhoSurf(i,j)
102     & + MAX( rhoMxL(i,j)*hMixCriteria, dRhoSmall )
103 jmc 1.3 hMixLayer(i,j,bi,bj) = rF(1)-R_low(i,j,bi,bj)
104 jmc 1.1 ENDDO
105     ENDDO
106     DO k = 2,Nr
107     C- potential density (reference level = surface level)
108 jmc 1.3 CALL FIND_RHO_2D(
109     I 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1,
110     I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
111     O rhoLoc,
112     I k, bi, bj, myThid )
113 jmc 1.1
114     DO j=1-Oly,sNy+Oly
115     DO i=1-Olx,sNx+Olx
116     IF ( k.LE.klowC(i,j,bi,bj) .AND.
117     & rhoLoc(i,j).GE.rhoMxL(i,j) ) THEN
118     IF ( rhoLoc(i,j).GT.rhoKm1(i,j) ) THEN
119     tmpFac = ( rhoMxL(i,j) - rhoKm1(i,j) )
120     & / ( rhoLoc(i,j) - rhoKm1(i,j) )
121     ELSE
122     tmpFac = 0.
123     ENDIF
124     hMixLayer(i,j,bi,bj) = rF(1)-rC(k-1)+tmpFac*drC(k)
125     rhoMxL(i,j) = rhoBigNb
126     ELSE
127     rhoKm1(i,j) = rhoLoc(i,j)
128     ENDIF
129     ENDDO
130     ENDDO
131     ENDDO
132    
133     ELSEIF ( method.EQ.2 ) THEN
134    
135     C-- Second method :
136     C where the local stratification exceed the mean stratification above
137     C (from surface down to here) by factor hMixCriteria
138    
139     c hMixCriteria = 2. _d 0
140 jmc 1.3 C-Note: dRhoSmall is hard coded for now but should become run-time parameter
141 jmc 1.2 dRhoSmall = 1. _d -2
142 jmc 1.1 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 jmc 1.2 DO k = 2,Nr-1
154 jmc 1.1 C- potential density (reference level = surface level)
155 jmc 1.3 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 jmc 1.1
161     DO j=1-Oly,sNy+Oly
162     DO i=1-Olx,sNx+Olx
163 jmc 1.2 IF ( k.LT.klowC(i,j,bi,bj)
164 jmc 1.1 & .AND. rhoMxL(i,j).GE.0. ) THEN
165 jmc 1.2 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 jmc 1.1 tmpFac = 0. _d 0
169     IF ( sigmAv.GT.0. _d 0 ) THEN
170 jmc 1.2 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 jmc 1.1 ENDIF
177 jmc 1.2 hMixLayer(i,j,bi,bj) = rF(1)-rF(k+1)
178     & - drF(k)*tmpFac*0.5 _d 0
179 jmc 1.1 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     #ifdef ALLOW_DIAGNOSTICS
191     IF ( useDiagnostics ) THEN
192     CALL DIAGNOSTICS_FILL( hMixLayer, 'MXLDEPTH',
193     & 0, 1, 1, bi, bj, myThid )
194     ENDIF
195     #endif /* ALLOW_DIAGNOSTICS */
196    
197     C-- end if calcMixLayerDepth
198     ENDIF
199    
200     RETURN
201     END

  ViewVC Help
Powered by ViewVC 1.1.22