1 |
C $Header: /u/gcmpack/MITgcm/pkg/my82/my82_calc.F,v 1.4 2005/05/30 07:41:45 mlosch 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, myTime, 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 |
c#include "DYNVARS.h" |
35 |
#include "MY82.h" |
36 |
c#include "FFIELDS.h" |
37 |
#include "GRID.h" |
38 |
|
39 |
C !INPUT PARAMETERS: =================================================== |
40 |
c Routine arguments |
41 |
c bi, bj - array indices on which to apply calculations |
42 |
c myTime - Current time in simulation |
43 |
|
44 |
INTEGER bi, bj |
45 |
_RL myTime |
46 |
INTEGER myThid |
47 |
|
48 |
#ifdef ALLOW_MY82 |
49 |
|
50 |
C !LOCAL VARIABLES: ==================================================== |
51 |
c Local constants |
52 |
C imin, imax, jmin, jmax - array computation indices |
53 |
C RiNumber - Richardson Number = -GH/GM |
54 |
C GH - buoyancy freqency after call to Ri_number, later |
55 |
C GH of M. Satoh, Eq. (11.3.45) |
56 |
C GM - vertical shear of velocity after call Ri_number, |
57 |
C later GM of M. Satoh, Eq. (11.3.44) |
58 |
INTEGER I, J, K |
59 |
INTEGER iMin ,iMax ,jMin ,jMax |
60 |
_RL RiFlux, RiTmp |
61 |
_RL SHtmp, bTmp, tkesquare, tkel |
62 |
_RL RiNumber(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
63 |
_RL GH(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
64 |
_RL GM(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
65 |
_RL SH(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
66 |
_RL SM(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
67 |
_RL tke(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
68 |
CEOP |
69 |
iMin = 2-OLx |
70 |
iMax = sNx+OLx-1 |
71 |
jMin = 2-OLy |
72 |
jMax = sNy+OLy-1 |
73 |
|
74 |
C Initialize local fields |
75 |
DO J=1-Oly,sNy+Oly |
76 |
DO I=1-Olx,sNx+Olx |
77 |
GH(I,J) = 0. _d 0 |
78 |
GM(I,J) = 0. _d 0 |
79 |
ENDDO |
80 |
ENDDO |
81 |
DO K = 1, Nr |
82 |
DO J=1-Oly,sNy+Oly |
83 |
DO I=1-Olx,sNx+Olx |
84 |
SH(I,J,K) = 0. _d 0 |
85 |
SM(I,J,K) = 0. _d 0 |
86 |
tke(I,J,K) = 0. _d 0 |
87 |
ENDDO |
88 |
ENDDO |
89 |
ENDDO |
90 |
C first k-loop |
91 |
C compute turbulent kinetic energy from richardson number |
92 |
DO K = 2, Nr |
93 |
CALL MY82_RI_NUMBER( |
94 |
I bi, bj, K, iMin, iMax, jMin, jMax, |
95 |
O RiNumber, GH, GM, |
96 |
I myTime, myThid ) |
97 |
DO J=jMin,jMax |
98 |
DO I=iMin,iMax |
99 |
RiTmp = MIN(RiNumber(I,J),RiMax) |
100 |
btmp = beta1+beta4*RiTmp |
101 |
C M. Satoh, Atmospheric Circulation Dynamics and General |
102 |
C Circulation models, Springer, 2004: Eq. (11.3.60) |
103 |
RiFlux =( btmp - SQRT(btmp*btmp - 4. _d 0 *beta2*beta3*RiTmp) ) |
104 |
& /(2. _d 0*beta2) |
105 |
C M. Satoh: Eq. (11.3.58) |
106 |
SHtmp = (alpha1-alpha2*RiFlux)/(1. _d 0-RiFlux) |
107 |
SH(I,J,K) = SHtmp |
108 |
SM(I,J,K) = SHtmp*(beta1-beta2*RiFlux)/(beta3-beta4*RiFlux) |
109 |
C M. Satoh: Eq. (11.3.53/55) |
110 |
tkesquare = MAX( 0. _d 0, |
111 |
& b1*(SH(I,J,K)*GH(I,J) + SM(I,J,K)*GM(I,J)) ) |
112 |
tke(I,J,K) = sqrt(tkesquare) |
113 |
CML if ( k.eq.2 .and. i.ge.1 .and. i.le.sNx .and. j.eq.1) |
114 |
CML & print '(A,3I3,8E13.5)', 'ml-my82', I,J,K, RiNumber(I,J), |
115 |
CML & RiTmp,RiFlux, |
116 |
CML & SH(I,J,K), SM(I,J,K), GH(I,J), GM(I,J), tke(I,J,K) |
117 |
ENDDO |
118 |
ENDDO |
119 |
C end of first k-loop |
120 |
ENDDO |
121 |
|
122 |
C re-initilialize GM and GH for abuse, they no longer have |
123 |
C the meaning of shear and negative buoyancy frequency |
124 |
DO J=jMin,jMax |
125 |
DO I=iMin,iMax |
126 |
GH(I,J) = 0. _d 0 |
127 |
GM(I,J) = 0. _d 0 |
128 |
ENDDO |
129 |
ENDDO |
130 |
C Find boundary length scale from energy weighted mean. |
131 |
C This is something like "center of mass" of the vertical |
132 |
C tke-distribution |
133 |
C begin second k-loop |
134 |
DO K = 2, Nr |
135 |
DO J=jMin,jMax |
136 |
DO I=iMin,iMax |
137 |
GM(I,J) = GM(I,J) + tke(I,J,K)*rF(K) |
138 |
GH(I,J) = GH(I,J) + tke(I,J,K) |
139 |
ENDDO |
140 |
ENDDO |
141 |
C end of second k-loop |
142 |
END DO |
143 |
C compute boundary length scale MYhbl |
144 |
DO J=jMin,jMax |
145 |
DO I=iMin,iMax |
146 |
IF ( GH(I,J) .EQ. 0. _d 0 ) THEN |
147 |
MYhbl(I,J,bi,bj) = 0. _d 0 |
148 |
ELSE |
149 |
MYhbl(I,J,bi,bj) = -GM(I,J)/GH(I,J)*MYhblScale |
150 |
ENDIF |
151 |
ENDDO |
152 |
ENDDO |
153 |
C begin third k-loop |
154 |
DO K = 1, Nr |
155 |
C integrate tke to find integral length scale |
156 |
DO J=jMin,jMax |
157 |
DO I=iMin,iMax |
158 |
tkel = MYhbl(I,J,bi,bj)*tke(I,J,K) |
159 |
C M. Satoh: Eq. (11.3.43) |
160 |
MYviscAr(I,J,K,bi,bj) = MYhbl(I,J,bi,bj)*tkel*SM(I,J,K) |
161 |
MYdiffKr(I,J,K,bi,bj) = MYhbl(I,J,bi,bj)*tkel*SH(I,J,K) |
162 |
C Set a minium (= background) value |
163 |
MYviscAr(I,J,K,bi,bj) = MAX(MYviscAr(I,J,K,bi,bj),viscAr) |
164 |
MYdiffKr(I,J,K,bi,bj) = MAX(MYdiffKr(I,J,K,bi,bj), |
165 |
& diffKrNrT(k) ) |
166 |
C Set a maximum and mask land point |
167 |
MYviscAr(I,J,K,bi,bj) = MIN(MYviscAr(I,J,K,bi,bj),MYviscMax) |
168 |
& * maskC(I,J,K,bi,bj) |
169 |
MYdiffKr(I,J,K,bi,bj) = MIN(MYdiffKr(I,J,K,bi,bj),MYdiffMax) |
170 |
& * maskC(I,J,K,bi,bj) |
171 |
ENDDO |
172 |
ENDDO |
173 |
C end third k-loop |
174 |
ENDDO |
175 |
|
176 |
#endif /* ALLOW_MY82 */ |
177 |
|
178 |
RETURN |
179 |
END |