/[MITgcm]/MITgcm_contrib/submesoscale/code/calc_oce_mxlayer.F
ViewVC logotype

Annotation of /MITgcm_contrib/submesoscale/code/calc_oce_mxlayer.F

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


Revision 1.2 - (hide annotations) (download)
Fri May 30 23:18:08 2008 UTC (17 years, 1 month ago) by dimitri
Branch: MAIN
Changes since 1.1: +1 -0 lines
Replacing KPPhbl with hMixLayer for submesoscale parameterization.

1 dimitri 1.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 dimitri 1.2 IF ( useGMRedi .AND. GM_SM_Ce.GT.0.0 ) calcMixLayerDepth = .TRUE.
71 dimitri 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     hMixLayer(i,j,bi,bj) = rF(1)-R_low(I,J,bi,bj)
104     ENDDO
105     ENDDO
106     DO k = 2,Nr
107     C- potential density (reference level = surface level)
108     CALL FIND_RHO(
109     I bi, bj, 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, K, 1,
110     I theta, salt,
111     O rhoLoc, 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(
155     I bi, bj, 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, K, 1,
156     I theta, salt,
157     O rhoLoc, myThid )
158    
159     DO j=1-Oly,sNy+Oly
160     DO i=1-Olx,sNx+Olx
161     IF ( k.LT.klowC(i,j,bi,bj)
162     & .AND. rhoMxL(i,j).GE.0. ) THEN
163     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     tmpFac = 0. _d 0
167     IF ( sigmAv.GT.0. _d 0 ) THEN
168     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     ENDIF
175     hMixLayer(i,j,bi,bj) = rF(1)-rF(k+1)
176     & - drF(k)*tmpFac*0.5 _d 0
177     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