/[MITgcm]/MITgcm/pkg/my82/my82_calc.F
ViewVC logotype

Contents of /MITgcm/pkg/my82/my82_calc.F

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


Revision 1.2 - (show annotations) (download)
Sun Oct 17 23:12:31 2004 UTC (19 years, 7 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint57d_post, checkpoint57b_post, checkpoint57c_pre, checkpoint55j_post, checkpoint56b_post, checkpoint57e_post, checkpoint55h_post, checkpoint56c_post, checkpoint57f_pre, checkpoint57a_post, checkpoint57a_pre, checkpoint55i_post, checkpoint57, checkpoint56, eckpoint57e_pre, checkpoint57c_post, checkpoint56a_post
Changes since 1.1: +3 -2 lines
allow to set a vertical profile of vertical diffusivity for T & S

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

  ViewVC Help
Powered by ViewVC 1.1.22