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