1 |
C $Header: /u/gcmpack/MITgcm/model/src/find_alpha.F,v 1.10 2002/09/18 16:30:34 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 |
|
14 |
C !DESCRIPTION: \bv |
15 |
C *==========================================================* |
16 |
C | o SUBROUTINE FIND_ALPHA |
17 |
C | Calculates [drho(S,T,z) / dT] of a horizontal slice |
18 |
C *==========================================================* |
19 |
C | |
20 |
C | k - is the Theta/Salt level |
21 |
C | kRef - determines pressure reference level |
22 |
C | (not used in 'LINEAR' mode) |
23 |
C | |
24 |
C | alphaloc - drho / dT (kg/m^3/C) |
25 |
C | |
26 |
C *==========================================================* |
27 |
C \ev |
28 |
|
29 |
C !USES: |
30 |
IMPLICIT NONE |
31 |
c Common |
32 |
#include "SIZE.h" |
33 |
#include "DYNVARS.h" |
34 |
#include "EEPARAMS.h" |
35 |
#include "PARAMS.h" |
36 |
#include "EOS.h" |
37 |
#include "GRID.h" |
38 |
|
39 |
C !INPUT/OUTPUT PARAMETERS: |
40 |
c Arguments |
41 |
integer bi,bj,iMin,iMax,jMin,jMax |
42 |
integer k ! Level of Theta/Salt slice |
43 |
integer kRef ! Pressure reference level |
44 |
_RL alphaloc(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
45 |
|
46 |
C !LOCAL VARIABLES: |
47 |
c Local |
48 |
integer i,j |
49 |
_RL refTemp,refSalt,tP,sP |
50 |
_RL t1, t2, t3, s1, s3o2, p1, p2, sp5, p1t1 |
51 |
_RL drhoP0dtheta, drhoP0dthetaFresh, drhoP0dthetaSalt |
52 |
_RL dKdtheta, dKdthetaFresh, dKdthetaSalt, dKdthetaPres |
53 |
_RL rhoP0(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
54 |
_RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
55 |
_RL dnum_dtheta, dden_dtheta |
56 |
_RL rhoDen(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
57 |
_RL rhoLoc(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
58 |
integer myThid |
59 |
CEOP |
60 |
|
61 |
Cml stop 'myThid is not properly defined in this subroutine' |
62 |
|
63 |
#ifdef CHECK_SALINITY_FOR_NEGATIVE_VALUES |
64 |
CALL LOOK_FOR_NEG_SALINITY( bi, bj, iMin, iMax, jMin, jMax, k, |
65 |
& sFld, myThid ) |
66 |
#endif |
67 |
|
68 |
if (equationOfState.eq.'LINEAR') then |
69 |
|
70 |
do j=jMin,jMax |
71 |
do i=iMin,iMax |
72 |
alphaloc(i,j) = -rhonil * tAlpha |
73 |
enddo |
74 |
enddo |
75 |
|
76 |
elseif (equationOfState.eq.'POLY3') then |
77 |
|
78 |
refTemp=eosRefT(kRef) |
79 |
refSalt=eosRefS(kRef) |
80 |
|
81 |
do j=jMin,jMax |
82 |
do i=iMin,iMax |
83 |
tP=theta(i,j,k,bi,bj)-refTemp |
84 |
sP=salt(i,j,k,bi,bj)-refSalt |
85 |
#ifdef USE_FACTORIZED_POLY |
86 |
alphaloc(i,j) = |
87 |
& ( eosC(6,kRef) |
88 |
& *tP*3. |
89 |
& +(eosC(7,kRef)*sP + eosC(3,kRef))*2. |
90 |
& )*tP |
91 |
& +(eosC(8,kRef)*sP + eosC(4,kRef) )*sP + eosC(1,kRef) |
92 |
& |
93 |
#else |
94 |
alphaloc(i,j) = |
95 |
& eosC(1,kRef) + |
96 |
& eosC(3,kRef)*tP*2. + |
97 |
& eosC(4,kRef) *sP + |
98 |
& eosC(6,kRef)*tP*tP*3. + |
99 |
& eosC(7,kRef)*tP*2. *sP + |
100 |
& eosC(8,kRef) *sP*sP |
101 |
#endif |
102 |
enddo |
103 |
enddo |
104 |
|
105 |
elseif ( equationOfState(1:5).eq.'JMD95' |
106 |
& .or. equationOfState.eq.'UNESCO' ) then |
107 |
C nonlinear equation of state in pressure coordinates |
108 |
|
109 |
CALL FIND_RHOP0( |
110 |
I bi, bj, iMin, iMax, jMin, jMax, k, |
111 |
I theta, salt, |
112 |
O rhoP0, |
113 |
I myThid ) |
114 |
|
115 |
CALL FIND_BULKMOD( |
116 |
I bi, bj, iMin, iMax, jMin, jMax, k, kRef, |
117 |
I theta, salt, |
118 |
O bulkMod, |
119 |
I myThid ) |
120 |
|
121 |
do j=jMin,jMax |
122 |
do i=iMin,iMax |
123 |
|
124 |
C abbreviations |
125 |
t1 = theta(i,j,k,bi,bj) |
126 |
t2 = t1*t1 |
127 |
t3 = t2*t1 |
128 |
|
129 |
s1 = salt(i,j,k,bi,bj) |
130 |
s3o2 = sqrt(s1*s1*s1) |
131 |
|
132 |
p1 = pressure(i,j,kRef,bi,bj)*SItoBar |
133 |
p2 = p1*p1 |
134 |
|
135 |
C d(rho)/d(theta) |
136 |
C of fresh water at p = 0 |
137 |
drhoP0dthetaFresh = |
138 |
& eosJMDCFw(2) |
139 |
& + 2.*eosJMDCFw(3)*t1 |
140 |
& + 3.*eosJMDCFw(4)*t2 |
141 |
& + 4.*eosJMDCFw(5)*t3 |
142 |
& + 5.*eosJMDCFw(6)*t3*t1 |
143 |
C of salt water at p = 0 |
144 |
drhoP0dthetaSalt = |
145 |
& s1*( |
146 |
& eosJMDCSw(2) |
147 |
& + 2.*eosJMDCSw(3)*t1 |
148 |
& + 3.*eosJMDCSw(4)*t2 |
149 |
& + 4.*eosJMDCSw(5)*t3 |
150 |
& ) |
151 |
& + s3o2*( |
152 |
& + eosJMDCSw(7) |
153 |
& + 2.*eosJMDCSw(8)*t1 |
154 |
& ) |
155 |
C d(bulk modulus)/d(theta) |
156 |
C of fresh water at p = 0 |
157 |
dKdthetaFresh = |
158 |
& eosJMDCKFw(2) |
159 |
& + 2.*eosJMDCKFw(3)*t1 |
160 |
& + 3.*eosJMDCKFw(4)*t2 |
161 |
& + 4.*eosJMDCKFw(5)*t3 |
162 |
C of sea water at p = 0 |
163 |
dKdthetaSalt = |
164 |
& s1*( eosJMDCKSw(2) |
165 |
& + 2.*eosJMDCKSw(3)*t1 |
166 |
& + 3.*eosJMDCKSw(4)*t2 |
167 |
& ) |
168 |
& + s3o2*( eosJMDCKSw(6) |
169 |
& + 2.*eosJMDCKSw(7)*t1 |
170 |
& ) |
171 |
C of sea water at p |
172 |
dKdthetaPres = |
173 |
& p1*( eosJMDCKP(2) |
174 |
& + 2.*eosJMDCKP(3)*t1 |
175 |
& + 3.*eosJMDCKP(4)*t2 |
176 |
& ) |
177 |
& + p1*s1*( eosJMDCKP(6) |
178 |
& + 2.*eosJMDCKP(7)*t1 |
179 |
& ) |
180 |
& + p2*( eosJMDCKP(10) |
181 |
& + 2.*eosJMDCKP(11)*t1 |
182 |
& ) |
183 |
& + p2*s1*( eosJMDCKP(13) |
184 |
& + 2.*eosJMDCKP(14)*t1 |
185 |
& ) |
186 |
|
187 |
drhoP0dtheta = drhoP0dthetaFresh |
188 |
& + drhoP0dthetaSalt |
189 |
dKdtheta = dKdthetaFresh |
190 |
& + dKdthetaSalt |
191 |
& + dKdthetaPres |
192 |
alphaloc(i,j) = |
193 |
& ( bulkmod(i,j)**2*drhoP0dtheta |
194 |
& - bulkmod(i,j)*p1*drhoP0dtheta |
195 |
& - rhoP0(i,j)*p1*dKdtheta ) |
196 |
& /( bulkmod(i,j) - p1 )**2 |
197 |
|
198 |
|
199 |
end do |
200 |
end do |
201 |
elseif ( equationOfState.eq.'MDJWF' ) then |
202 |
|
203 |
CALL FIND_RHO( bi, bj, iMin, iMax, jMin, jMax, k, kRef, |
204 |
& theta, salt, rhoLoc, myThid ) |
205 |
CALL FIND_RHODEN( bi, bj, iMin, iMax, jMin, jMax, k, kRef, |
206 |
& theta, salt, rhoDen, myThid ) |
207 |
|
208 |
do j=jMin,jMax |
209 |
do i=iMin,iMax |
210 |
t1 = theta(i,j,k,bi,bj) |
211 |
t2 = t1*t1 |
212 |
s1 = salt(i,j,k,bi,bj) |
213 |
sp5 = sqrt(s1) |
214 |
|
215 |
p1 = pressure(i,j,kRef,bi,bj)*SItodBar |
216 |
p1t1 = p1*t1 |
217 |
|
218 |
dnum_dtheta = eosMDJWFnum(1) |
219 |
& + t1*(2.*eosMDJWFnum(2) + 3.*eosMDJWFnum(3)*t1) |
220 |
& + eosMDJWFnum(5)*s1 |
221 |
& + p1t1*(2.*eosMDJWFnum(8) + 2.*eosMDJWFnum(11)*p1) |
222 |
|
223 |
dden_dtheta = eosMDJWFden(1) |
224 |
& + t1*(2.*eosMDJWFden(2) |
225 |
& + t1*(3.*eosMDJWFden(3) |
226 |
& + 4.*eosMDJWFden(4)*t1 ) ) |
227 |
& + s1*(eosMDJWFden(6) |
228 |
& + t1*(3.*eosMDJWFden(7)*t1 |
229 |
& + 2.*eosMDJWFden(9)*sp5 ) ) |
230 |
& + p1*p1*(3.*eosMDJWFden(11)*t2 + eosMDJWFden(12)*p1) |
231 |
|
232 |
alphaLoc(i,j) = rhoDen(i,j)*(dnum_dtheta |
233 |
& - rhoLoc(i,j)*dden_dtheta) |
234 |
|
235 |
end do |
236 |
end do |
237 |
|
238 |
else |
239 |
write(*,*) 'FIND_ALPHA: equationOfState = ',equationOfState |
240 |
stop 'FIND_ALPHA: "equationOfState" has illegal value' |
241 |
endif |
242 |
|
243 |
return |
244 |
end |
245 |
|
246 |
subroutine FIND_BETA ( |
247 |
I bi, bj, iMin, iMax, jMin, jMax, k, kRef, |
248 |
O betaloc ) |
249 |
C /==========================================================\ |
250 |
C | o SUBROUTINE FIND_BETA | |
251 |
C | Calculates [drho(S,T,z) / dS] of a horizontal slice | |
252 |
C |==========================================================| |
253 |
C | | |
254 |
C | k - is the Theta/Salt level | |
255 |
C | kRef - determines pressure reference level | |
256 |
C | (not used in 'LINEAR' mode) | |
257 |
C | | |
258 |
C | betaloc - drho / dS (kg/m^3/PSU) | |
259 |
C | | |
260 |
C \==========================================================/ |
261 |
implicit none |
262 |
|
263 |
c Common |
264 |
#include "SIZE.h" |
265 |
#include "DYNVARS.h" |
266 |
#include "EEPARAMS.h" |
267 |
#include "PARAMS.h" |
268 |
#include "EOS.h" |
269 |
#include "GRID.h" |
270 |
|
271 |
c Arguments |
272 |
integer bi,bj,iMin,iMax,jMin,jMax |
273 |
integer k ! Level of Theta/Salt slice |
274 |
integer kRef ! Pressure reference level |
275 |
_RL betaloc(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
276 |
|
277 |
c Local |
278 |
integer i,j |
279 |
_RL refTemp,refSalt,tP,sP |
280 |
_RL t1, t2, t3, s1, s3o2, p1, sp5, p1t1 |
281 |
_RL drhoP0dS |
282 |
_RL dKdS, dKdSSalt, dKdSPres |
283 |
_RL rhoP0(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
284 |
_RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
285 |
_RL dnum_dsalt, dden_dsalt |
286 |
_RL rhoDen(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
287 |
_RL rhoLoc(1-Olx:sNx+Olx,1-Oly:sNy+Oly) |
288 |
integer myThid |
289 |
CEOP |
290 |
|
291 |
Cml stop 'myThid is not properly defined in this subroutine' |
292 |
|
293 |
#ifdef CHECK_SALINITY_FOR_NEGATIVE_VALUES |
294 |
CALL LOOK_FOR_NEG_SALINITY( bi, bj, iMin, iMax, jMin, jMax, k, |
295 |
& sFld, myThid ) |
296 |
#endif |
297 |
|
298 |
if (equationOfState.eq.'LINEAR') then |
299 |
|
300 |
do j=jMin,jMax |
301 |
do i=iMin,iMax |
302 |
betaloc(i,j) = rhonil * sBeta |
303 |
enddo |
304 |
enddo |
305 |
|
306 |
elseif (equationOfState.eq.'POLY3') then |
307 |
|
308 |
refTemp=eosRefT(kRef) |
309 |
refSalt=eosRefS(kRef) |
310 |
|
311 |
do j=jMin,jMax |
312 |
do i=iMin,iMax |
313 |
tP=theta(i,j,k,bi,bj)-refTemp |
314 |
sP=salt(i,j,k,bi,bj)-refSalt |
315 |
#ifdef USE_FACTORIZED_POLY |
316 |
betaloc(i,j) = |
317 |
& ( eosC(9,kRef)*sP*3. + eosC(5,kRef)*2. )*sP + eosC(2,kRef) |
318 |
& + ( eosC(7,kRef)*tP |
319 |
& +eosC(8,kRef)*sP*2. + eosC(4,kRef) |
320 |
& )*tP |
321 |
#else |
322 |
betaloc(i,j) = |
323 |
& eosC(2,kRef) + |
324 |
& eosC(4,kRef)*tP + |
325 |
& eosC(5,kRef) *sP*2. + |
326 |
& eosC(7,kRef)*tP*tP + |
327 |
& eosC(8,kRef)*tP *sP*2. + |
328 |
& eosC(9,kRef) *sP*sP*3. |
329 |
#endif |
330 |
enddo |
331 |
enddo |
332 |
|
333 |
elseif ( equationOfState(1:5).eq.'JMD95' |
334 |
& .or. equationOfState.eq.'UNESCO' ) then |
335 |
C nonlinear equation of state in pressure coordinates |
336 |
|
337 |
CALL FIND_RHOP0( |
338 |
I bi, bj, iMin, iMax, jMin, jMax, k, |
339 |
I theta, salt, |
340 |
O rhoP0, |
341 |
I myThid ) |
342 |
|
343 |
CALL FIND_BULKMOD( |
344 |
I bi, bj, iMin, iMax, jMin, jMax, k, kRef, |
345 |
I theta, salt, |
346 |
O bulkMod, |
347 |
I myThid ) |
348 |
|
349 |
do j=jMin,jMax |
350 |
do i=iMin,iMax |
351 |
|
352 |
C abbreviations |
353 |
t1 = theta(i,j,k,bi,bj) |
354 |
t2 = t1*t1 |
355 |
t3 = t2*t1 |
356 |
|
357 |
s1 = salt(i,j,k,bi,bj) |
358 |
s3o2 = 1.5*sqrt(s1) |
359 |
|
360 |
p1 = pressure(i,j,kRef,bi,bj)*SItoBar |
361 |
|
362 |
C d(rho)/d(S) |
363 |
C of fresh water at p = 0 |
364 |
drhoP0dS = 0. _d 0 |
365 |
C of salt water at p = 0 |
366 |
drhoP0dS = drhoP0dS |
367 |
& + eosJMDCSw(1) |
368 |
& + eosJMDCSw(2)*t1 |
369 |
& + eosJMDCSw(3)*t2 |
370 |
& + eosJMDCSw(4)*t3 |
371 |
& + eosJMDCSw(5)*t3*t1 |
372 |
& + s3o2*( |
373 |
& eosJMDCSw(6) |
374 |
& + eosJMDCSw(7)*t1 |
375 |
& + eosJMDCSw(8)*t2 |
376 |
& ) |
377 |
& + 2*eosJMDCSw(9)*s1 |
378 |
C d(bulk modulus)/d(S) |
379 |
C of fresh water at p = 0 |
380 |
dKdS = 0. _d 0 |
381 |
C of sea water at p = 0 |
382 |
dKdSSalt = |
383 |
& eosJMDCKSw(1) |
384 |
& + eosJMDCKSw(2)*t1 |
385 |
& + eosJMDCKSw(3)*t2 |
386 |
& + eosJMDCKSw(4)*t3 |
387 |
& + s3o2*( eosJMDCKSw(5) |
388 |
& + eosJMDCKSw(6)*t1 |
389 |
& + eosJMDCKSw(7)*t2 |
390 |
& ) |
391 |
|
392 |
C of sea water at p |
393 |
dKdSPres = |
394 |
& p1*( eosJMDCKP(5) |
395 |
& + eosJMDCKP(6)*t1 |
396 |
& + eosJMDCKP(7)*t2 |
397 |
& ) |
398 |
& + s3o2*p1*eosJMDCKP(8) |
399 |
& + p1*p1*( eosJMDCKP(12) |
400 |
& + eosJMDCKP(13)*t1 |
401 |
& + eosJMDCKP(14)*t2 |
402 |
& ) |
403 |
|
404 |
dKdS = dKdSSalt + dKdSPres |
405 |
|
406 |
betaloc(i,j) = |
407 |
& ( bulkmod(i,j)**2*drhoP0dS |
408 |
& - bulkmod(i,j)*p1*drhoP0dS |
409 |
& - rhoP0(i,j)*p1*dKdS ) |
410 |
& /( bulkmod(i,j) - p1 )**2 |
411 |
|
412 |
|
413 |
end do |
414 |
end do |
415 |
elseif ( equationOfState.eq.'MDJWF' ) then |
416 |
|
417 |
CALL FIND_RHO( bi, bj, iMin, iMax, jMin, jMax, k, kRef, |
418 |
& theta, salt, rhoLoc, myThid ) |
419 |
CALL FIND_RHODEN( bi, bj, iMin, iMax, jMin, jMax, k, kRef, |
420 |
& theta, salt, rhoDen, myThid ) |
421 |
|
422 |
do j=jMin,jMax |
423 |
do i=iMin,iMax |
424 |
t1 = theta(i,j,k,bi,bj) |
425 |
t2 = t1*t1 |
426 |
s1 = salt(i,j,k,bi,bj) |
427 |
sp5 = sqrt(s1) |
428 |
|
429 |
p1 = pressure(i,j,k,bi,bj)*SItodBar |
430 |
p1t1 = p1*t1 |
431 |
|
432 |
dnum_dsalt = eosMDJWFnum(4) |
433 |
& + eosMDJWFnum(5)*t1 |
434 |
& + 2.*eosMDJWFnum(6)*s1 + eosMDJWFnum(9)*p1 |
435 |
dden_dsalt = eosMDJWFden(5) |
436 |
& + t1*( eosMDJWFden(6) + eosMDJWFden(7)*t2 ) |
437 |
& + 1.5*sp5*(eosMDJWFden(8) + eosMDJWFden(9)*t2) |
438 |
|
439 |
betaLoc(i,j) = rhoDen(i,j)*( dnum_dsalt |
440 |
& - rhoLoc(i,j)*dden_dsalt ) |
441 |
|
442 |
end do |
443 |
end do |
444 |
|
445 |
else |
446 |
write(*,*) 'FIND_BETA: equationOfState = ',equationOfState |
447 |
stop 'FIND_BETA: "equationOfState" has illegal value' |
448 |
endif |
449 |
|
450 |
return |
451 |
end |