1 |
C $Header: /u/gcmpack/MITgcm/model/src/find_alpha.F,v 1.20 2012/08/12 20:24:23 jmc Exp $ |
2 |
C $Name: $ |
3 |
|
4 |
#include "PACKAGES_CONFIG.h" |
5 |
#include "CPP_OPTIONS.h" |
6 |
#define USE_FACTORIZED_POLY |
7 |
|
8 |
C-- File find_alpha.F: Calculate expansion Coeff |
9 |
C-- Contents |
10 |
C-- o FIND_ALPHA |
11 |
C-- o FIND_BETA |
12 |
|
13 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
14 |
CBOP |
15 |
C !ROUTINE: FIND_ALPHA |
16 |
C !INTERFACE: |
17 |
SUBROUTINE FIND_ALPHA ( |
18 |
I bi, bj, iMin, iMax, jMin, jMax, k, kRef, |
19 |
O alphaLoc, |
20 |
I myThid ) |
21 |
|
22 |
C !DESCRIPTION: \bv |
23 |
C *==========================================================* |
24 |
C | o SUBROUTINE FIND_ALPHA |
25 |
C | Calculates [drho(S,T,z) / dT] of a horizontal slice |
26 |
C *==========================================================* |
27 |
C | k - is the Theta/Salt level |
28 |
C | kRef - determines pressure reference level |
29 |
C | (not used in 'LINEAR' mode) |
30 |
C | alphaLoc - drho / dT (kg/m^3/C) |
31 |
C *==========================================================* |
32 |
C \ev |
33 |
|
34 |
C !USES: |
35 |
IMPLICIT NONE |
36 |
C === Global variables === |
37 |
#include "SIZE.h" |
38 |
#include "DYNVARS.h" |
39 |
#include "EEPARAMS.h" |
40 |
#include "PARAMS.h" |
41 |
#include "EOS.h" |
42 |
#include "GRID.h" |
43 |
|
44 |
C !INPUT/OUTPUT PARAMETERS: |
45 |
C == Routine arguments == |
46 |
C k :: Level of Theta/Salt slice |
47 |
C kRef :: Pressure reference level |
48 |
c myThid :: thread number for this instance of the routine |
49 |
INTEGER myThid |
50 |
INTEGER bi,bj,iMin,iMax,jMin,jMax |
51 |
INTEGER k |
52 |
INTEGER kRef |
53 |
_RL alphaLoc(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
54 |
|
55 |
C !LOCAL VARIABLES: |
56 |
C == Local variables == |
57 |
INTEGER i,j |
58 |
_RL refTemp,refSalt,tP,sP |
59 |
_RL t1, t2, t3, s1, s3o2, p1, p2, sp5, p1t1 |
60 |
_RL ct, sa, sqrtsa, p |
61 |
_RL drhoP0dtheta, drhoP0dthetaFresh, drhoP0dthetaSalt |
62 |
_RL dKdtheta, dKdthetaFresh, dKdthetaSalt, dKdthetaPres |
63 |
_RL locPres(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
64 |
_RL rhoP0 (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
65 |
_RL bulkMod(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
66 |
_RL dnum_dtheta, dden_dtheta |
67 |
_RL rhoDen (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
68 |
_RL rhoLoc (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
69 |
CEOP |
70 |
|
71 |
#ifdef CHECK_SALINITY_FOR_NEGATIVE_VALUES |
72 |
c CALL LOOK_FOR_NEG_SALINITY( |
73 |
c I iMin, iMax, jMin, jMax, |
74 |
c U sFld, |
75 |
c I k, bi, bj, myThid ) |
76 |
#endif |
77 |
|
78 |
IF (equationOfState.EQ.'LINEAR') THEN |
79 |
|
80 |
DO j=jMin,jMax |
81 |
DO i=iMin,iMax |
82 |
alphaLoc(i,j) = -rhonil * tAlpha |
83 |
ENDDO |
84 |
ENDDO |
85 |
|
86 |
ELSEIF (equationOfState.EQ.'POLY3') THEN |
87 |
|
88 |
refTemp=eosRefT(kRef) |
89 |
refSalt=eosRefS(kRef) |
90 |
|
91 |
DO j=jMin,jMax |
92 |
DO i=iMin,iMax |
93 |
tP=theta(i,j,k,bi,bj)-refTemp |
94 |
sP=salt(i,j,k,bi,bj)-refSalt |
95 |
#ifdef USE_FACTORIZED_POLY |
96 |
alphaLoc(i,j) = |
97 |
& ( eosC(6,kRef) |
98 |
& *tP*3. |
99 |
& +(eosC(7,kRef)*sP + eosC(3,kRef))*2. |
100 |
& )*tP |
101 |
& +(eosC(8,kRef)*sP + eosC(4,kRef) )*sP + eosC(1,kRef) |
102 |
& |
103 |
#else |
104 |
alphaLoc(i,j) = |
105 |
& eosC(1,kRef) + |
106 |
& eosC(3,kRef)*tP*2. + |
107 |
& eosC(4,kRef) *sP + |
108 |
& eosC(6,kRef)*tP*tP*3. + |
109 |
& eosC(7,kRef)*tP*2. *sP + |
110 |
& eosC(8,kRef) *sP*sP |
111 |
#endif |
112 |
ENDDO |
113 |
ENDDO |
114 |
|
115 |
ELSEIF ( equationOfState(1:5).EQ.'JMD95' |
116 |
& .OR. equationOfState.EQ.'UNESCO' ) THEN |
117 |
C nonlinear equation of state in pressure coordinates |
118 |
|
119 |
CALL PRESSURE_FOR_EOS( |
120 |
I bi, bj, iMin, iMax, jMin, jMax, kRef, |
121 |
O locPres, |
122 |
I myThid ) |
123 |
|
124 |
CALL FIND_RHOP0( |
125 |
I iMin, iMax, jMin, jMax, |
126 |
I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj), |
127 |
O rhoP0, |
128 |
I myThid ) |
129 |
|
130 |
CALL FIND_BULKMOD( |
131 |
I iMin, iMax, jMin, jMax, locPres, |
132 |
I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj), |
133 |
O bulkMod, |
134 |
I myThid ) |
135 |
|
136 |
DO j=jMin,jMax |
137 |
DO i=iMin,iMax |
138 |
|
139 |
C abbreviations |
140 |
t1 = theta(i,j,k,bi,bj) |
141 |
t2 = t1*t1 |
142 |
t3 = t2*t1 |
143 |
|
144 |
#if (defined ALLOW_AUTODIFF && defined TARGET_NEC_SX) |
145 |
C this line makes TAF create vectorizable code |
146 |
s3o2 = 0. _d 0 |
147 |
#endif |
148 |
s1 = salt(i,j,k,bi,bj) |
149 |
IF ( s1 .GT. 0. _d 0 ) THEN |
150 |
s3o2 = SQRT(s1*s1*s1) |
151 |
ELSE |
152 |
s1 = 0. _d 0 |
153 |
s3o2 = 0. _d 0 |
154 |
ENDIF |
155 |
|
156 |
p1 = locPres(i,j)*SItoBar |
157 |
p2 = p1*p1 |
158 |
|
159 |
C d(rho)/d(theta) |
160 |
C of fresh water at p = 0 |
161 |
drhoP0dthetaFresh = |
162 |
& eosJMDCFw(2) |
163 |
& + 2.*eosJMDCFw(3)*t1 |
164 |
& + 3.*eosJMDCFw(4)*t2 |
165 |
& + 4.*eosJMDCFw(5)*t3 |
166 |
& + 5.*eosJMDCFw(6)*t3*t1 |
167 |
C of salt water at p = 0 |
168 |
drhoP0dthetaSalt = |
169 |
& s1*( |
170 |
& eosJMDCSw(2) |
171 |
& + 2.*eosJMDCSw(3)*t1 |
172 |
& + 3.*eosJMDCSw(4)*t2 |
173 |
& + 4.*eosJMDCSw(5)*t3 |
174 |
& ) |
175 |
& + s3o2*( |
176 |
& + eosJMDCSw(7) |
177 |
& + 2.*eosJMDCSw(8)*t1 |
178 |
& ) |
179 |
C d(bulk modulus)/d(theta) |
180 |
C of fresh water at p = 0 |
181 |
dKdthetaFresh = |
182 |
& eosJMDCKFw(2) |
183 |
& + 2.*eosJMDCKFw(3)*t1 |
184 |
& + 3.*eosJMDCKFw(4)*t2 |
185 |
& + 4.*eosJMDCKFw(5)*t3 |
186 |
C of sea water at p = 0 |
187 |
dKdthetaSalt = |
188 |
& s1*( eosJMDCKSw(2) |
189 |
& + 2.*eosJMDCKSw(3)*t1 |
190 |
& + 3.*eosJMDCKSw(4)*t2 |
191 |
& ) |
192 |
& + s3o2*( eosJMDCKSw(6) |
193 |
& + 2.*eosJMDCKSw(7)*t1 |
194 |
& ) |
195 |
C of sea water at p |
196 |
dKdthetaPres = |
197 |
& p1*( eosJMDCKP(2) |
198 |
& + 2.*eosJMDCKP(3)*t1 |
199 |
& + 3.*eosJMDCKP(4)*t2 |
200 |
& ) |
201 |
& + p1*s1*( eosJMDCKP(6) |
202 |
& + 2.*eosJMDCKP(7)*t1 |
203 |
& ) |
204 |
& + p2*( eosJMDCKP(10) |
205 |
& + 2.*eosJMDCKP(11)*t1 |
206 |
& ) |
207 |
& + p2*s1*( eosJMDCKP(13) |
208 |
& + 2.*eosJMDCKP(14)*t1 |
209 |
& ) |
210 |
|
211 |
drhoP0dtheta = drhoP0dthetaFresh |
212 |
& + drhoP0dthetaSalt |
213 |
dKdtheta = dKdthetaFresh |
214 |
& + dKdthetaSalt |
215 |
& + dKdthetaPres |
216 |
alphaLoc(i,j) = |
217 |
& ( bulkmod(i,j)**2*drhoP0dtheta |
218 |
& - bulkmod(i,j)*p1*drhoP0dtheta |
219 |
& - rhoP0(i,j)*p1*dKdtheta ) |
220 |
& /( bulkmod(i,j) - p1 )**2 |
221 |
|
222 |
ENDDO |
223 |
ENDDO |
224 |
ELSEIF ( equationOfState.EQ.'MDJWF' ) THEN |
225 |
|
226 |
CALL PRESSURE_FOR_EOS( |
227 |
I bi, bj, iMin, iMax, jMin, jMax, kRef, |
228 |
O locPres, |
229 |
I myThid ) |
230 |
|
231 |
CALL FIND_RHONUM( |
232 |
I iMin, iMax, jMin, jMax, locPres, |
233 |
I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj), |
234 |
O rhoLoc, |
235 |
I myThid ) |
236 |
|
237 |
CALL FIND_RHODEN( |
238 |
I iMin, iMax, jMin, jMax, locPres, |
239 |
I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj), |
240 |
O rhoDen, |
241 |
I myThid ) |
242 |
|
243 |
DO j=jMin,jMax |
244 |
DO i=iMin,iMax |
245 |
t1 = theta(i,j,k,bi,bj) |
246 |
t2 = t1*t1 |
247 |
s1 = salt(i,j,k,bi,bj) |
248 |
#if (defined ALLOW_AUTODIFF && defined TARGET_NEC_SX) |
249 |
C this line makes TAF create vectorizable code |
250 |
sp5 = 0. _d 0 |
251 |
#endif |
252 |
IF ( s1 .GT. 0. _d 0 ) THEN |
253 |
sp5 = SQRT(s1) |
254 |
ELSE |
255 |
s1 = 0. _d 0 |
256 |
sp5 = 0. _d 0 |
257 |
ENDIF |
258 |
|
259 |
p1 = locPres(i,j)*SItodBar |
260 |
p1t1 = p1*t1 |
261 |
|
262 |
dnum_dtheta = eosMDJWFnum(1) |
263 |
& + t1*(2.*eosMDJWFnum(2) + 3.*eosMDJWFnum(3)*t1) |
264 |
& + eosMDJWFnum(5)*s1 |
265 |
& + p1t1*(2.*eosMDJWFnum(8) + 2.*eosMDJWFnum(11)*p1) |
266 |
|
267 |
dden_dtheta = eosMDJWFden(1) |
268 |
& + t1*(2.*eosMDJWFden(2) |
269 |
& + t1*(3.*eosMDJWFden(3) |
270 |
& + 4.*eosMDJWFden(4)*t1 ) ) |
271 |
& + s1*(eosMDJWFden(6) |
272 |
& + t1*(3.*eosMDJWFden(7)*t1 |
273 |
& + 2.*eosMDJWFden(9)*sp5 ) ) |
274 |
& + p1*p1*(3.*eosMDJWFden(11)*t2 + eosMDJWFden(12)*p1) |
275 |
|
276 |
alphaLoc(i,j) = rhoDen(i,j)*(dnum_dtheta |
277 |
& - (rhoLoc(i,j)*rhoDen(i,j))*dden_dtheta) |
278 |
|
279 |
ENDDO |
280 |
ENDDO |
281 |
|
282 |
ELSEIF ( equationOfState.EQ.'TEOS10' ) THEN |
283 |
|
284 |
CALL PRESSURE_FOR_EOS( |
285 |
I bi, bj, iMin, iMax, jMin, jMax, kRef, |
286 |
O locPres, |
287 |
I myThid ) |
288 |
|
289 |
CALL FIND_RHOTEOS( |
290 |
I iMin, iMax, jMin, jMax, locPres, |
291 |
I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj), |
292 |
O rhoLoc, rhoDen, |
293 |
I myThid ) |
294 |
|
295 |
DO j=jMin,jMax |
296 |
DO i=iMin,iMax |
297 |
ct = theta(i,j,k,bi,bj) |
298 |
sa = salt(i,j,k,bi,bj) |
299 |
#if (defined ALLOW_AUTODIFF && defined TARGET_NEC_SX) |
300 |
C this line makes TAF create vectorizable code |
301 |
sqrtsa = 0. _d 0 |
302 |
#endif |
303 |
IF ( sa .GT. 0. _d 0 ) THEN |
304 |
sqrtsa = SQRT(sa) |
305 |
ELSE |
306 |
sa = 0. _d 0 |
307 |
sqrtsa = 0. _d 0 |
308 |
ENDIF |
309 |
p = locPres(i,j)*SItodBar |
310 |
|
311 |
dnum_dtheta = teos(02) |
312 |
& + ct*(2.*teos(03) + 3.*teos(04)*ct) |
313 |
& + sa*(teos(06) + 2.*teos(07)*ct |
314 |
& + sqrtsa*(teos(09) + ct*(2.*teos(10) + 3.*teos(11)*ct))) |
315 |
& + p*( teos(13) + 2.*teos(14)*ct + sa*2.*teos(16) |
316 |
& + p*(teos(18) + 2.*teos(19)*ct)) |
317 |
|
318 |
dden_dtheta = teos(22) |
319 |
& + ct*(2.*teos(23) + ct*(3.*teos(24) + 4.*teos(25)*ct)) |
320 |
& + sa*(teos(27) |
321 |
& + ct*(2.*teos(28) + ct*(3.*teos(29) + 4.*teos(30)*ct)) |
322 |
& + sqrtsa*(teos(32) |
323 |
& + ct*(2.*teos(33) + ct*(3.*teos(34) + 4.*teos(35)*ct)))) |
324 |
& + p*(teos(38) + ct*(2.*teos(39) + 3.*teos(40)*ct) |
325 |
& + teos(42) |
326 |
& + p*(teos(44) + 2.*teos(45)*ct + teos(46)*sa |
327 |
& + p*teos(48) )) |
328 |
|
329 |
alphaLoc(i,j) = rhoDen(i,j)*(dnum_dtheta |
330 |
& - (rhoLoc(i,j)*rhoDen(i,j))*dden_dtheta) |
331 |
|
332 |
ENDDO |
333 |
ENDDO |
334 |
|
335 |
ELSE |
336 |
WRITE(*,*) 'FIND_ALPHA: equationOfState = ',equationOfState |
337 |
STOP 'FIND_ALPHA: "equationOfState" has illegal value' |
338 |
ENDIF |
339 |
|
340 |
RETURN |
341 |
END |
342 |
|
343 |
C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----| |
344 |
CBOP |
345 |
C !ROUTINE: FIND_BETA |
346 |
C !INTERFACE: |
347 |
SUBROUTINE FIND_BETA ( |
348 |
I bi, bj, iMin, iMax, jMin, jMax, k, kRef, |
349 |
O betaLoc, |
350 |
I myThid ) |
351 |
|
352 |
C !DESCRIPTION: \bv |
353 |
C *==========================================================* |
354 |
C | o SUBROUTINE FIND_BETA |
355 |
C | Calculates [drho(S,T,z) / dS] of a horizontal slice |
356 |
C |==========================================================* |
357 |
C | k - is the Theta/Salt level |
358 |
C | kRef - determines pressure reference level |
359 |
C | (not used in 'LINEAR' mode) |
360 |
C | betaLoc - drho / dS (kg/m^3/PSU) |
361 |
C *==========================================================* |
362 |
C \ev |
363 |
|
364 |
C !USES: |
365 |
IMPLICIT NONE |
366 |
C === Global variables === |
367 |
#include "SIZE.h" |
368 |
#include "DYNVARS.h" |
369 |
#include "EEPARAMS.h" |
370 |
#include "PARAMS.h" |
371 |
#include "EOS.h" |
372 |
#include "GRID.h" |
373 |
|
374 |
C !INPUT/OUTPUT PARAMETERS: |
375 |
C == Routine arguments == |
376 |
C k :: Level of Theta/Salt slice |
377 |
C kRef :: Pressure reference level |
378 |
c myThid :: thread number for this instance of the routine |
379 |
INTEGER myThid |
380 |
INTEGER bi,bj,iMin,iMax,jMin,jMax |
381 |
INTEGER k |
382 |
INTEGER kRef |
383 |
_RL betaLoc(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
384 |
|
385 |
C !LOCAL VARIABLES: |
386 |
C == Local variables == |
387 |
INTEGER i,j |
388 |
_RL refTemp,refSalt,tP,sP |
389 |
_RL t1, t2, t3, s1, s3o2, p1, sp5, p1t1 |
390 |
_RL ct, sa, sqrtsa, p |
391 |
_RL drhoP0dS |
392 |
_RL dKdS, dKdSSalt, dKdSPres |
393 |
_RL locPres(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
394 |
_RL rhoP0 (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
395 |
_RL bulkMod(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
396 |
_RL dnum_dsalt, dden_dsalt |
397 |
_RL rhoDen (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
398 |
_RL rhoLoc (1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
399 |
CEOP |
400 |
|
401 |
#ifdef CHECK_SALINITY_FOR_NEGATIVE_VALUES |
402 |
c CALL LOOK_FOR_NEG_SALINITY( |
403 |
c I iMin, iMax, jMin, jMax, |
404 |
c U sFld, |
405 |
c I k, bi, bj, myThid ) |
406 |
#endif |
407 |
|
408 |
IF (equationOfState.EQ.'LINEAR') THEN |
409 |
|
410 |
DO j=jMin,jMax |
411 |
DO i=iMin,iMax |
412 |
betaLoc(i,j) = rhonil * sBeta |
413 |
ENDDO |
414 |
ENDDO |
415 |
|
416 |
ELSEIF (equationOfState.EQ.'POLY3') THEN |
417 |
|
418 |
refTemp=eosRefT(kRef) |
419 |
refSalt=eosRefS(kRef) |
420 |
|
421 |
DO j=jMin,jMax |
422 |
DO i=iMin,iMax |
423 |
tP=theta(i,j,k,bi,bj)-refTemp |
424 |
sP=salt(i,j,k,bi,bj)-refSalt |
425 |
#ifdef USE_FACTORIZED_POLY |
426 |
betaLoc(i,j) = |
427 |
& ( eosC(9,kRef)*sP*3. + eosC(5,kRef)*2. )*sP + eosC(2,kRef) |
428 |
& + ( eosC(7,kRef)*tP |
429 |
& +eosC(8,kRef)*sP*2. + eosC(4,kRef) |
430 |
& )*tP |
431 |
#else |
432 |
betaLoc(i,j) = |
433 |
& eosC(2,kRef) + |
434 |
& eosC(4,kRef)*tP + |
435 |
& eosC(5,kRef) *sP*2. + |
436 |
& eosC(7,kRef)*tP*tP + |
437 |
& eosC(8,kRef)*tP *sP*2. + |
438 |
& eosC(9,kRef) *sP*sP*3. |
439 |
#endif |
440 |
ENDDO |
441 |
ENDDO |
442 |
|
443 |
ELSEIF ( equationOfState(1:5).EQ.'JMD95' |
444 |
& .OR. equationOfState.EQ.'UNESCO' ) THEN |
445 |
C nonlinear equation of state in pressure coordinates |
446 |
|
447 |
CALL PRESSURE_FOR_EOS( |
448 |
I bi, bj, iMin, iMax, jMin, jMax, kRef, |
449 |
O locPres, |
450 |
I myThid ) |
451 |
|
452 |
CALL FIND_RHOP0( |
453 |
I iMin, iMax, jMin, jMax, |
454 |
I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj), |
455 |
O rhoP0, |
456 |
I myThid ) |
457 |
|
458 |
CALL FIND_BULKMOD( |
459 |
I iMin, iMax, jMin, jMax, locPres, |
460 |
I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj), |
461 |
O bulkMod, |
462 |
I myThid ) |
463 |
|
464 |
DO j=jMin,jMax |
465 |
DO i=iMin,iMax |
466 |
|
467 |
C abbreviations |
468 |
t1 = theta(i,j,k,bi,bj) |
469 |
t2 = t1*t1 |
470 |
t3 = t2*t1 |
471 |
|
472 |
s1 = salt(i,j,k,bi,bj) |
473 |
#if (defined ALLOW_AUTODIFF && defined TARGET_NEC_SX) |
474 |
C this line makes TAF create vectorizable code |
475 |
s3o2 = 0. _d 0 |
476 |
#endif |
477 |
IF ( s1 .GT. 0. _d 0 ) THEN |
478 |
s3o2 = 1.5*SQRT(s1) |
479 |
ELSE |
480 |
s1 = 0. _d 0 |
481 |
s3o2 = 0. _d 0 |
482 |
ENDIF |
483 |
|
484 |
p1 = locPres(i,j)*SItoBar |
485 |
|
486 |
C d(rho)/d(S) |
487 |
C of fresh water at p = 0 |
488 |
drhoP0dS = 0. _d 0 |
489 |
C of salt water at p = 0 |
490 |
drhoP0dS = drhoP0dS |
491 |
& + eosJMDCSw(1) |
492 |
& + eosJMDCSw(2)*t1 |
493 |
& + eosJMDCSw(3)*t2 |
494 |
& + eosJMDCSw(4)*t3 |
495 |
& + eosJMDCSw(5)*t3*t1 |
496 |
& + s3o2*( |
497 |
& eosJMDCSw(6) |
498 |
& + eosJMDCSw(7)*t1 |
499 |
& + eosJMDCSw(8)*t2 |
500 |
& ) |
501 |
& + 2*eosJMDCSw(9)*s1 |
502 |
C d(bulk modulus)/d(S) |
503 |
C of fresh water at p = 0 |
504 |
dKdS = 0. _d 0 |
505 |
C of sea water at p = 0 |
506 |
dKdSSalt = |
507 |
& eosJMDCKSw(1) |
508 |
& + eosJMDCKSw(2)*t1 |
509 |
& + eosJMDCKSw(3)*t2 |
510 |
& + eosJMDCKSw(4)*t3 |
511 |
& + s3o2*( eosJMDCKSw(5) |
512 |
& + eosJMDCKSw(6)*t1 |
513 |
& + eosJMDCKSw(7)*t2 |
514 |
& ) |
515 |
|
516 |
C of sea water at p |
517 |
dKdSPres = |
518 |
& p1*( eosJMDCKP(5) |
519 |
& + eosJMDCKP(6)*t1 |
520 |
& + eosJMDCKP(7)*t2 |
521 |
& ) |
522 |
& + s3o2*p1*eosJMDCKP(8) |
523 |
& + p1*p1*( eosJMDCKP(12) |
524 |
& + eosJMDCKP(13)*t1 |
525 |
& + eosJMDCKP(14)*t2 |
526 |
& ) |
527 |
|
528 |
dKdS = dKdSSalt + dKdSPres |
529 |
|
530 |
betaLoc(i,j) = |
531 |
& ( bulkmod(i,j)**2*drhoP0dS |
532 |
& - bulkmod(i,j)*p1*drhoP0dS |
533 |
& - rhoP0(i,j)*p1*dKdS ) |
534 |
& /( bulkmod(i,j) - p1 )**2 |
535 |
|
536 |
ENDDO |
537 |
ENDDO |
538 |
ELSEIF ( equationOfState.EQ.'MDJWF' ) THEN |
539 |
|
540 |
CALL PRESSURE_FOR_EOS( |
541 |
I bi, bj, iMin, iMax, jMin, jMax, kRef, |
542 |
O locPres, |
543 |
I myThid ) |
544 |
|
545 |
CALL FIND_RHONUM( |
546 |
I iMin, iMax, jMin, jMax, locPres, |
547 |
I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj), |
548 |
O rhoLoc, |
549 |
I myThid ) |
550 |
|
551 |
CALL FIND_RHODEN( |
552 |
I iMin, iMax, jMin, jMax, locPres, |
553 |
I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj), |
554 |
O rhoDen, |
555 |
I myThid ) |
556 |
|
557 |
DO j=jMin,jMax |
558 |
DO i=iMin,iMax |
559 |
t1 = theta(i,j,k,bi,bj) |
560 |
t2 = t1*t1 |
561 |
s1 = salt(i,j,k,bi,bj) |
562 |
#if (defined ALLOW_AUTODIFF && defined TARGET_NEC_SX) |
563 |
C this line makes TAF create vectorizable code |
564 |
sp5 = 0. _d 0 |
565 |
#endif |
566 |
IF ( s1 .GT. 0. _d 0 ) THEN |
567 |
sp5 = SQRT(s1) |
568 |
ELSE |
569 |
s1 = 0. _d 0 |
570 |
sp5 = 0. _d 0 |
571 |
ENDIF |
572 |
|
573 |
p1 = locPres(i,j)*SItodBar |
574 |
p1t1 = p1*t1 |
575 |
|
576 |
dnum_dsalt = eosMDJWFnum(4) |
577 |
& + eosMDJWFnum(5)*t1 |
578 |
& + 2.*eosMDJWFnum(6)*s1 + eosMDJWFnum(9)*p1 |
579 |
dden_dsalt = eosMDJWFden(5) |
580 |
& + t1*( eosMDJWFden(6) + eosMDJWFden(7)*t2 ) |
581 |
& + 1.5*sp5*(eosMDJWFden(8) + eosMDJWFden(9)*t2) |
582 |
|
583 |
betaLoc(i,j) = rhoDen(i,j)*( dnum_dsalt |
584 |
& - (rhoLoc(i,j)*rhoDen(i,j))*dden_dsalt ) |
585 |
|
586 |
ENDDO |
587 |
ENDDO |
588 |
|
589 |
ELSEIF ( equationOfState.EQ.'TEOS10' ) THEN |
590 |
|
591 |
CALL PRESSURE_FOR_EOS( |
592 |
I bi, bj, iMin, iMax, jMin, jMax, kRef, |
593 |
O locPres, |
594 |
I myThid ) |
595 |
|
596 |
CALL FIND_RHOTEOS( |
597 |
I iMin, iMax, jMin, jMax, locPres, |
598 |
I theta(1-OLx,1-OLy,k,bi,bj), salt(1-OLx,1-OLy,k,bi,bj), |
599 |
O rhoLoc, rhoDen, |
600 |
I myThid ) |
601 |
|
602 |
DO j=jMin,jMax |
603 |
DO i=iMin,iMax |
604 |
ct = theta(i,j,k,bi,bj) |
605 |
sa = salt(i,j,k,bi,bj) |
606 |
#if (defined ALLOW_AUTODIFF && defined TARGET_NEC_SX) |
607 |
C this line makes TAF create vectorizable code |
608 |
sqrtsa = 0. _d 0 |
609 |
#endif |
610 |
IF ( sa .GT. 0. _d 0 ) THEN |
611 |
sqrtsa = SQRT(sa) |
612 |
ELSE |
613 |
sa = 0. _d 0 |
614 |
sqrtsa = 0. _d 0 |
615 |
ENDIF |
616 |
p = locPres(i,j)*SItodBar |
617 |
|
618 |
dnum_dsalt = teos(05) + ct*(teos(06) + teos(07)*ct) |
619 |
& + 1.5*sqrtsa*(teos(08) |
620 |
& + ct*(teos(09) + ct*(teos(10) + teos(11)*ct))) |
621 |
& + p*(teos(15) + teos(16)*ct + p*teos(20)) |
622 |
|
623 |
dden_dsalt = teos(26) |
624 |
& + ct*(teos(27) + ct*(teos(28) + ct*(teos(29) + teos(30)*ct))) |
625 |
& + 2.*teos(36)*sa |
626 |
& + 1.5*sqrtsa*(teos(31) + ct*(teos(32) + ct*(teos(33) |
627 |
& + ct*(teos(34) + teos(35)*ct)))) |
628 |
& + p*(teos(41) + teos(42)*ct + p*teos(46)) |
629 |
|
630 |
betaLoc(i,j) = rhoDen(i,j)*( dnum_dsalt |
631 |
& - (rhoLoc(i,j)*rhoDen(i,j))*dden_dsalt ) |
632 |
|
633 |
ENDDO |
634 |
ENDDO |
635 |
|
636 |
ELSE |
637 |
WRITE(*,*) 'FIND_BETA: equationOfState = ',equationOfState |
638 |
STOP 'FIND_BETA: "equationOfState" has illegal value' |
639 |
ENDIF |
640 |
|
641 |
RETURN |
642 |
END |