1 |
C $Header: /u/gcmpack/MITgcm/pkg/my82/my82_calc.F,v 1.6 2009/10/08 20:07:53 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "MY82_OPTIONS.h" |
5 |
|
6 |
CBOP |
7 |
C !ROUTINE: MY82_CALC |
8 |
|
9 |
C !INTERFACE: ====================================================== |
10 |
subroutine MY82_CALC( |
11 |
I bi, bj, sigmaR, myTime, myIter, myThid ) |
12 |
|
13 |
C !DESCRIPTION: \bv |
14 |
C *==========================================================* |
15 |
C | SUBROUTINE MY82_CALC | |
16 |
C | o Compute all MY82 fields defined in MY82.h | |
17 |
C *==========================================================* |
18 |
C | This subroutine is based on SPEM code | |
19 |
C *==========================================================* |
20 |
IMPLICIT NONE |
21 |
C |
22 |
C-------------------------------------------------------------------- |
23 |
|
24 |
C global parameters updated by pp_calc |
25 |
C PPviscAz - PP eddy viscosity coefficient (m^2/s) |
26 |
C PPdiffKzT - PP diffusion coefficient for temperature (m^2/s) |
27 |
C |
28 |
C \ev |
29 |
|
30 |
C !USES: ============================================================ |
31 |
#include "SIZE.h" |
32 |
#include "EEPARAMS.h" |
33 |
#include "PARAMS.h" |
34 |
#include "MY82.h" |
35 |
#include "GRID.h" |
36 |
|
37 |
C !INPUT PARAMETERS: =================================================== |
38 |
C Routine arguments |
39 |
C bi, bj :: Current tile indices |
40 |
C sigmaR :: Vertical gradient of iso-neutral density |
41 |
C myTime :: Current time in simulation |
42 |
C myIter :: Current time-step number |
43 |
C myThid :: My Thread Id number |
44 |
INTEGER bi, bj |
45 |
_RL sigmaR(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
46 |
_RL myTime |
47 |
INTEGER myIter |
48 |
INTEGER myThid |
49 |
|
50 |
#ifdef ALLOW_MY82 |
51 |
C !LOCAL VARIABLES: ==================================================== |
52 |
c Local constants |
53 |
C imin, imax, jmin, jmax - array computation indices |
54 |
C RiNumber - Richardson Number = -GH/GM |
55 |
C GH - buoyancy freqency after call to Ri_number, later |
56 |
C GH of M. Satoh, Eq. (11.3.45) |
57 |
C GM - vertical shear of velocity after call Ri_number, |
58 |
C later GM of M. Satoh, Eq. (11.3.44) |
59 |
INTEGER I, J, K |
60 |
INTEGER iMin ,iMax ,jMin ,jMax |
61 |
_RL RiFlux, RiTmp |
62 |
_RL SHtmp, bTmp, tkesquare, tkel |
63 |
_RL RiNumber(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
64 |
_RL GH(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
65 |
_RL GM(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
66 |
_RL SH(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
67 |
_RL SM(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
68 |
_RL tke(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
69 |
CEOP |
70 |
iMin = 2-OLx |
71 |
iMax = sNx+OLx-1 |
72 |
jMin = 2-OLy |
73 |
jMax = sNy+OLy-1 |
74 |
|
75 |
C Initialize local fields |
76 |
DO J=1-OLy,sNy+OLy |
77 |
DO I=1-OLx,sNx+OLx |
78 |
GH(I,J) = 0. _d 0 |
79 |
GM(I,J) = 0. _d 0 |
80 |
ENDDO |
81 |
ENDDO |
82 |
DO K = 1, Nr |
83 |
DO J=1-OLy,sNy+OLy |
84 |
DO I=1-OLx,sNx+OLx |
85 |
SH(I,J,K) = 0. _d 0 |
86 |
SM(I,J,K) = 0. _d 0 |
87 |
tke(I,J,K) = 0. _d 0 |
88 |
ENDDO |
89 |
ENDDO |
90 |
ENDDO |
91 |
C first k-loop |
92 |
C compute turbulent kinetic energy from richardson number |
93 |
DO K = 2, Nr |
94 |
CALL MY82_RI_NUMBER( |
95 |
I bi, bj, K, iMin, iMax, jMin, jMax, |
96 |
O RiNumber, GH, GM, |
97 |
I myTime, myThid ) |
98 |
DO J=jMin,jMax |
99 |
DO I=iMin,iMax |
100 |
RiTmp = MIN(RiNumber(I,J),RiMax) |
101 |
btmp = beta1+beta4*RiTmp |
102 |
C M. Satoh, Atmospheric Circulation Dynamics and General |
103 |
C Circulation models, Springer, 2004: Eq. (11.3.60) |
104 |
RiFlux =( btmp - SQRT(btmp*btmp - 4. _d 0 *beta2*beta3*RiTmp) ) |
105 |
& /(2. _d 0*beta2) |
106 |
C M. Satoh: Eq. (11.3.58) |
107 |
SHtmp = (alpha1-alpha2*RiFlux)/(1. _d 0-RiFlux) |
108 |
SH(I,J,K) = SHtmp |
109 |
SM(I,J,K) = SHtmp*(beta1-beta2*RiFlux)/(beta3-beta4*RiFlux) |
110 |
C M. Satoh: Eq. (11.3.53/55) |
111 |
tkesquare = MAX( 0. _d 0, |
112 |
& b1*(SH(I,J,K)*GH(I,J) + SM(I,J,K)*GM(I,J)) ) |
113 |
tke(I,J,K) = sqrt(tkesquare) |
114 |
CML if ( k.eq.2 .and. i.ge.1 .and. i.le.sNx .and. j.eq.1) |
115 |
CML & print '(A,3I3,8E13.5)', 'ml-my82', I,J,K, RiNumber(I,J), |
116 |
CML & RiTmp,RiFlux, |
117 |
CML & SH(I,J,K), SM(I,J,K), GH(I,J), GM(I,J), tke(I,J,K) |
118 |
ENDDO |
119 |
ENDDO |
120 |
C end of first k-loop |
121 |
ENDDO |
122 |
|
123 |
C re-initilialize GM and GH for abuse, they no longer have |
124 |
C the meaning of shear and negative buoyancy frequency |
125 |
DO J=jMin,jMax |
126 |
DO I=iMin,iMax |
127 |
GH(I,J) = 0. _d 0 |
128 |
GM(I,J) = 0. _d 0 |
129 |
ENDDO |
130 |
ENDDO |
131 |
C Find boundary length scale from energy weighted mean. |
132 |
C This is something like "center of mass" of the vertical |
133 |
C tke-distribution |
134 |
C begin second k-loop |
135 |
DO K = 2, Nr |
136 |
DO J=jMin,jMax |
137 |
DO I=iMin,iMax |
138 |
GM(I,J) = GM(I,J) + tke(I,J,K)*rF(K) |
139 |
GH(I,J) = GH(I,J) + tke(I,J,K) |
140 |
ENDDO |
141 |
ENDDO |
142 |
C end of second k-loop |
143 |
END DO |
144 |
C compute boundary length scale MYhbl |
145 |
DO J=jMin,jMax |
146 |
DO I=iMin,iMax |
147 |
IF ( GH(I,J) .EQ. 0. _d 0 ) THEN |
148 |
MYhbl(I,J,bi,bj) = 0. _d 0 |
149 |
ELSE |
150 |
MYhbl(I,J,bi,bj) = -GM(I,J)/GH(I,J)*MYhblScale |
151 |
ENDIF |
152 |
ENDDO |
153 |
ENDDO |
154 |
C begin third k-loop |
155 |
DO K = 1, Nr |
156 |
C integrate tke to find integral length scale |
157 |
DO J=jMin,jMax |
158 |
DO I=iMin,iMax |
159 |
tkel = MYhbl(I,J,bi,bj)*tke(I,J,K) |
160 |
C M. Satoh: Eq. (11.3.43) |
161 |
MYviscAr(I,J,K,bi,bj) = MYhbl(I,J,bi,bj)*tkel*SM(I,J,K) |
162 |
MYdiffKr(I,J,K,bi,bj) = MYhbl(I,J,bi,bj)*tkel*SH(I,J,K) |
163 |
C Set a minium (= background) value |
164 |
MYviscAr(I,J,K,bi,bj) = MAX(MYviscAr(I,J,K,bi,bj), |
165 |
& viscArnr(k) ) |
166 |
MYdiffKr(I,J,K,bi,bj) = MAX(MYdiffKr(I,J,K,bi,bj), |
167 |
& diffKrNrT(k) ) |
168 |
C Set a maximum and mask land point |
169 |
MYviscAr(I,J,K,bi,bj) = MIN(MYviscAr(I,J,K,bi,bj),MYviscMax) |
170 |
& * maskC(I,J,K,bi,bj) |
171 |
MYdiffKr(I,J,K,bi,bj) = MIN(MYdiffKr(I,J,K,bi,bj),MYdiffMax) |
172 |
& * maskC(I,J,K,bi,bj) |
173 |
ENDDO |
174 |
ENDDO |
175 |
C end third k-loop |
176 |
ENDDO |
177 |
|
178 |
#endif /* ALLOW_MY82 */ |
179 |
|
180 |
RETURN |
181 |
END |