/[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.1 - (hide annotations) (download)
Sun Jun 3 22:50:47 2007 UTC (16 years, 11 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint59c
add subroutine to calculate ocean mixed-layer depth.

1 jmc 1.1 C $Header: /u/gcmpack/MITgcm/model/src/calc_grad_phi_surf.F,v 1.2 2001/09/26 18:09:13 cnh 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(
108     I bi, bj, 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, K, 1,
109     I theta, salt,
110     O rhoLoc, myThid )
111    
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     DO j=1-Oly,sNy+Oly
139     DO i=1-Olx,sNx+Olx
140     IF ( klowC(i,j,bi,bj) .GT. 0 ) THEN
141     hMixLayer(i,j,bi,bj) = drF(1)
142     rhoMxL(i,j) = 1.
143     ELSE
144     hMixLayer(i,j,bi,bj) = rF(1)
145     rhoMxL(i,j) = -1.
146     ENDIF
147     ENDDO
148     ENDDO
149     DO k = 3,Nr
150     C- potential density (reference level = surface level)
151     CALL FIND_RHO(
152     I bi, bj, 1-OLx, sNx+OLx, 1-OLy, sNy+OLy, K, 1,
153     I theta, salt,
154     O rhoLoc, myThid )
155    
156     DO j=1-Oly,sNy+Oly
157     DO i=1-Olx,sNx+Olx
158     IF ( k.LE.klowC(i,j,bi,bj)
159     & .AND. rhoMxL(i,j).GE.0. ) THEN
160     sigmAv = ( rhoLoc(i,j)-rhoSurf(i,j) ) / ( rC(1)-rC(k) )
161     IF ( -sigmaR(i,j,k).GE.sigmAv*hMixCriteria ) THEN
162     tmpFac = 0. _d 0
163     IF ( sigmAv.GT.0. _d 0 ) THEN
164     tmpFac = sigmaR(i,j,k-1)/sigmaR(i,j,k)
165     tmpFac = ( tmpFac*hMixCriteria - 1. _d 0 )
166     & / ( hMixCriteria - 1. _d 0 )
167     tmpFac = MAX( 0. _d 0, MIN( tmpFac, 1. _d 0 ) )
168     ENDIF
169     hMixLayer(i,j,bi,bj) = rF(1)-rC(k-1)+tmpFac*drC(k)
170     rhoMxL(i,j) = -1.
171     ENDIF
172     ENDIF
173     ENDDO
174     ENDDO
175     ENDDO
176    
177     ELSE
178     STOP 'S/R CALC_OCE_MXLAYER: invalid method'
179     ENDIF
180    
181     #ifdef ALLOW_DIAGNOSTICS
182     IF ( useDiagnostics ) THEN
183     CALL DIAGNOSTICS_FILL( hMixLayer, 'MXLDEPTH',
184     & 0, 1, 1, bi, bj, myThid )
185     ENDIF
186     #endif /* ALLOW_DIAGNOSTICS */
187    
188     C-- end if calcMixLayerDepth
189     ENDIF
190    
191     RETURN
192     END

  ViewVC Help
Powered by ViewVC 1.1.22