/[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.4 - (hide annotations) (download)
Tue Aug 12 22:26:39 2008 UTC (15 years, 9 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint61d, checkpoint61e, checkpoint61c
Changes since 1.3: +1 -2 lines
SURFACE.h no longer needed.

1 jmc 1.4 C $Header: /u/gcmpack/MITgcm/model/src/calc_oce_mxlayer.F,v 1.3 2008/08/11 22:25:52 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     IF ( useGMRedi ) calcMixLayerDepth = GM_taper_scheme.EQ.'fm07'
69     #endif
70     #ifdef ALLOW_DIAGNOSTICS
71     IF ( useDiagnostics.AND. .NOT.calcMixLayerDepth ) THEN
72     calcMixLayerDepth = DIAGNOSTICS_IS_ON('MXLDEPTH',myThid)
73     ENDIF
74     #endif
75     IF ( calcMixLayerDepth ) THEN
76    
77     C-- Select which "method" to use:
78     method = 0
79     IF ( hMixCriteria.LT.0. ) method = 1
80     IF ( hMixCriteria.GT.1. ) method = 2
81    
82     IF ( method.EQ.1 ) THEN
83    
84     C-- First method :
85     C where the potential density (ref.lev=surface) is larger than
86     C surface density plus Delta_Rho = hMixCriteria * Alpha(surf)
87     C = density of water which is -hMixCriteria colder than surface water
88    
89     c hMixCriteria = -0.8 _d 0
90     dRhoSmall = 1. _d -6
91     rhoBigNb = rhoConst*1. _d 10
92     CALL FIND_ALPHA(
93     I bi, bj, 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1, 1,
94     O rhoMxL, myThid )
95    
96     DO j=1-Oly,sNy+Oly
97     DO i=1-Olx,sNx+Olx
98     rhoKm1(i,j) = rhoSurf(i,j)
99     rhoMxL(i,j) = rhoSurf(i,j)
100     & + MAX( rhoMxL(i,j)*hMixCriteria, dRhoSmall )
101 jmc 1.3 hMixLayer(i,j,bi,bj) = rF(1)-R_low(i,j,bi,bj)
102 jmc 1.1 ENDDO
103     ENDDO
104     DO k = 2,Nr
105     C- potential density (reference level = surface level)
106 jmc 1.3 CALL FIND_RHO_2D(
107     I 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1,
108     I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
109     O rhoLoc,
110     I k, bi, bj, myThid )
111 jmc 1.1
112     DO j=1-Oly,sNy+Oly
113     DO i=1-Olx,sNx+Olx
114     IF ( k.LE.klowC(i,j,bi,bj) .AND.
115     & rhoLoc(i,j).GE.rhoMxL(i,j) ) THEN
116     IF ( rhoLoc(i,j).GT.rhoKm1(i,j) ) THEN
117     tmpFac = ( rhoMxL(i,j) - rhoKm1(i,j) )
118     & / ( rhoLoc(i,j) - rhoKm1(i,j) )
119     ELSE
120     tmpFac = 0.
121     ENDIF
122     hMixLayer(i,j,bi,bj) = rF(1)-rC(k-1)+tmpFac*drC(k)
123     rhoMxL(i,j) = rhoBigNb
124     ELSE
125     rhoKm1(i,j) = rhoLoc(i,j)
126     ENDIF
127     ENDDO
128     ENDDO
129     ENDDO
130    
131     ELSEIF ( method.EQ.2 ) THEN
132    
133     C-- Second method :
134     C where the local stratification exceed the mean stratification above
135     C (from surface down to here) by factor hMixCriteria
136    
137     c hMixCriteria = 2. _d 0
138 jmc 1.3 C-Note: dRhoSmall is hard coded for now but should become run-time parameter
139 jmc 1.2 dRhoSmall = 1. _d -2
140 jmc 1.1 DO j=1-Oly,sNy+Oly
141     DO i=1-Olx,sNx+Olx
142     IF ( klowC(i,j,bi,bj) .GT. 0 ) THEN
143     hMixLayer(i,j,bi,bj) = drF(1)
144     rhoMxL(i,j) = 1.
145     ELSE
146     hMixLayer(i,j,bi,bj) = rF(1)
147     rhoMxL(i,j) = -1.
148     ENDIF
149     ENDDO
150     ENDDO
151 jmc 1.2 DO k = 2,Nr-1
152 jmc 1.1 C- potential density (reference level = surface level)
153 jmc 1.3 CALL FIND_RHO_2D(
154     I 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, 1,
155     I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj),
156     O rhoLoc,
157     I k, bi, bj, myThid )
158 jmc 1.1
159     DO j=1-Oly,sNy+Oly
160     DO i=1-Olx,sNx+Olx
161 jmc 1.2 IF ( k.LT.klowC(i,j,bi,bj)
162 jmc 1.1 & .AND. rhoMxL(i,j).GE.0. ) THEN
163 jmc 1.2 sigmAv = ( rhoLoc(i,j)-rhoSurf(i,j)+dRhoSmall )
164     & / ( rC(1)-rC(k) )
165     IF ( -sigmaR(i,j,k+1).GT.sigmAv*hMixCriteria ) THEN
166 jmc 1.1 tmpFac = 0. _d 0
167     IF ( sigmAv.GT.0. _d 0 ) THEN
168 jmc 1.2 tmpFac = hMixCriteria*sigmaR(i,j,k)/sigmaR(i,j,k+1)
169     IF ( tmpFac .GT. 1. _d 0 ) THEN
170     tmpFac = 1. _d 0
171     & + ( tmpFac - 1. _d 0 )/( hMixCriteria - 1. _d 0 )
172     ENDIF
173     tmpFac = MAX( 0. _d 0, MIN( tmpFac, 2. _d 0 ) )
174 jmc 1.1 ENDIF
175 jmc 1.2 hMixLayer(i,j,bi,bj) = rF(1)-rF(k+1)
176     & - drF(k)*tmpFac*0.5 _d 0
177 jmc 1.1 rhoMxL(i,j) = -1.
178     ENDIF
179     ENDIF
180     ENDDO
181     ENDDO
182     ENDDO
183    
184     ELSE
185     STOP 'S/R CALC_OCE_MXLAYER: invalid method'
186     ENDIF
187    
188     #ifdef ALLOW_DIAGNOSTICS
189     IF ( useDiagnostics ) THEN
190     CALL DIAGNOSTICS_FILL( hMixLayer, 'MXLDEPTH',
191     & 0, 1, 1, bi, bj, myThid )
192     ENDIF
193     #endif /* ALLOW_DIAGNOSTICS */
194    
195     C-- end if calcMixLayerDepth
196     ENDIF
197    
198     RETURN
199     END

  ViewVC Help
Powered by ViewVC 1.1.22