/[MITgcm]/MITgcm/model/src/find_alpha.F
ViewVC logotype

Contents of /MITgcm/model/src/find_alpha.F

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


Revision 1.16 - (show annotations) (download)
Tue May 1 04:09:36 2007 UTC (17 years, 1 month ago) by mlosch
Branch: MAIN
CVS Tags: checkpoint60, checkpoint61, checkpoint59q, checkpoint59p, checkpoint59r, checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint59c, checkpoint59b, checkpoint59m, checkpoint59l, checkpoint59o, checkpoint59n, checkpoint59i, checkpoint59h, checkpoint59k, checkpoint59j, checkpoint61b, checkpoint61a
Changes since 1.15: +13 -9 lines
 add myThid to all kpp routines (long overdue)

1 C $Header: /u/gcmpack/MITgcm/model/src/find_alpha.F,v 1.15 2003/10/27 22:32:55 heimbach 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 drhoP0dtheta, drhoP0dthetaFresh, drhoP0dthetaSalt
57 _RL dKdtheta, dKdthetaFresh, dKdthetaSalt, dKdthetaPres
58 _RL locPres(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
59 _RL rhoP0 (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
60 _RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
61 _RL dnum_dtheta, dden_dtheta
62 _RL rhoDen (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
63 _RL rhoLoc (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
64 CEOP
65
66 #ifdef CHECK_SALINITY_FOR_NEGATIVE_VALUES
67 CALL LOOK_FOR_NEG_SALINITY( bi, bj, iMin, iMax, jMin, jMax, k,
68 & sFld, myThid )
69 #endif
70
71 IF (equationOfState.EQ.'LINEAR') THEN
72
73 DO j=jMin,jMax
74 DO i=iMin,iMax
75 alphaloc(i,j) = -rhonil * tAlpha
76 ENDDO
77 ENDDO
78
79 ELSEIF (equationOfState.EQ.'POLY3') THEN
80
81 refTemp=eosRefT(kRef)
82 refSalt=eosRefS(kRef)
83
84 DO j=jMin,jMax
85 DO i=iMin,iMax
86 tP=theta(i,j,k,bi,bj)-refTemp
87 sP=salt(i,j,k,bi,bj)-refSalt
88 #ifdef USE_FACTORIZED_POLY
89 alphaloc(i,j) =
90 & ( eosC(6,kRef)
91 & *tP*3.
92 & +(eosC(7,kRef)*sP + eosC(3,kRef))*2.
93 & )*tP
94 & +(eosC(8,kRef)*sP + eosC(4,kRef) )*sP + eosC(1,kRef)
95 &
96 #else
97 alphaloc(i,j) =
98 & eosC(1,kRef) +
99 & eosC(3,kRef)*tP*2. +
100 & eosC(4,kRef) *sP +
101 & eosC(6,kRef)*tP*tP*3. +
102 & eosC(7,kRef)*tP*2. *sP +
103 & eosC(8,kRef) *sP*sP
104 #endif
105 ENDDO
106 ENDDO
107
108 ELSEIF ( equationOfState(1:5).EQ.'JMD95'
109 & .OR. equationOfState.EQ.'UNESCO' ) THEN
110 C nonlinear equation of state in pressure coordinates
111
112 CALL PRESSURE_FOR_EOS(
113 I bi, bj, iMin, iMax, jMin, jMax, kRef,
114 O locPres,
115 I myThid )
116
117 CALL FIND_RHOP0(
118 I bi, bj, iMin, iMax, jMin, jMax, k,
119 I theta, salt,
120 O rhoP0,
121 I myThid )
122
123 CALL FIND_BULKMOD(
124 I bi, bj, iMin, iMax, jMin, jMax, k,
125 I locPres, theta, salt,
126 O bulkMod,
127 I myThid )
128
129 DO j=jMin,jMax
130 DO i=iMin,iMax
131
132 C abbreviations
133 t1 = theta(i,j,k,bi,bj)
134 t2 = t1*t1
135 t3 = t2*t1
136
137 s1 = salt(i,j,k,bi,bj)
138 IF ( s1 .GT. 0. _d 0 ) THEN
139 s3o2 = SQRT(s1*s1*s1)
140 ELSE
141 s1 = 0. _d 0
142 s3o2 = 0. _d 0
143 ENDIF
144
145 p1 = locPres(i,j)*SItoBar
146 p2 = p1*p1
147
148 C d(rho)/d(theta)
149 C of fresh water at p = 0
150 drhoP0dthetaFresh =
151 & eosJMDCFw(2)
152 & + 2.*eosJMDCFw(3)*t1
153 & + 3.*eosJMDCFw(4)*t2
154 & + 4.*eosJMDCFw(5)*t3
155 & + 5.*eosJMDCFw(6)*t3*t1
156 C of salt water at p = 0
157 drhoP0dthetaSalt =
158 & s1*(
159 & eosJMDCSw(2)
160 & + 2.*eosJMDCSw(3)*t1
161 & + 3.*eosJMDCSw(4)*t2
162 & + 4.*eosJMDCSw(5)*t3
163 & )
164 & + s3o2*(
165 & + eosJMDCSw(7)
166 & + 2.*eosJMDCSw(8)*t1
167 & )
168 C d(bulk modulus)/d(theta)
169 C of fresh water at p = 0
170 dKdthetaFresh =
171 & eosJMDCKFw(2)
172 & + 2.*eosJMDCKFw(3)*t1
173 & + 3.*eosJMDCKFw(4)*t2
174 & + 4.*eosJMDCKFw(5)*t3
175 C of sea water at p = 0
176 dKdthetaSalt =
177 & s1*( eosJMDCKSw(2)
178 & + 2.*eosJMDCKSw(3)*t1
179 & + 3.*eosJMDCKSw(4)*t2
180 & )
181 & + s3o2*( eosJMDCKSw(6)
182 & + 2.*eosJMDCKSw(7)*t1
183 & )
184 C of sea water at p
185 dKdthetaPres =
186 & p1*( eosJMDCKP(2)
187 & + 2.*eosJMDCKP(3)*t1
188 & + 3.*eosJMDCKP(4)*t2
189 & )
190 & + p1*s1*( eosJMDCKP(6)
191 & + 2.*eosJMDCKP(7)*t1
192 & )
193 & + p2*( eosJMDCKP(10)
194 & + 2.*eosJMDCKP(11)*t1
195 & )
196 & + p2*s1*( eosJMDCKP(13)
197 & + 2.*eosJMDCKP(14)*t1
198 & )
199
200 drhoP0dtheta = drhoP0dthetaFresh
201 & + drhoP0dthetaSalt
202 dKdtheta = dKdthetaFresh
203 & + dKdthetaSalt
204 & + dKdthetaPres
205 alphaloc(i,j) =
206 & ( bulkmod(i,j)**2*drhoP0dtheta
207 & - bulkmod(i,j)*p1*drhoP0dtheta
208 & - rhoP0(i,j)*p1*dKdtheta )
209 & /( bulkmod(i,j) - p1 )**2
210
211
212 ENDDO
213 ENDDO
214 ELSEIF ( equationOfState.EQ.'MDJWF' ) THEN
215
216 CALL PRESSURE_FOR_EOS(
217 I bi, bj, iMin, iMax, jMin, jMax, kRef,
218 O locPres,
219 I myThid )
220
221 c$$$ CALL FIND_RHO( bi, bj, iMin, iMax, jMin, jMax, k,
222 c$$$ & kRef, theta, salt, rhoLoc, myThid )
223
224 CALL FIND_RHONUM( bi, bj, iMin, iMax, jMin, jMax, k,
225 & locPres, theta, salt, rhoLoc, myThid )
226
227 CALL FIND_RHODEN( bi, bj, iMin, iMax, jMin, jMax, k,
228 & locPres, theta, salt, rhoDen, myThid )
229
230 DO j=jMin,jMax
231 DO i=iMin,iMax
232 t1 = theta(i,j,k,bi,bj)
233 t2 = t1*t1
234 s1 = salt(i,j,k,bi,bj)
235 IF ( s1 .GT. 0. _d 0 ) THEN
236 sp5 = SQRT(s1)
237 ELSE
238 s1 = 0. _d 0
239 sp5 = 0. _d 0
240 ENDIF
241
242 p1 = locPres(i,j)*SItodBar
243 p1t1 = p1*t1
244
245 dnum_dtheta = eosMDJWFnum(1)
246 & + t1*(2.*eosMDJWFnum(2) + 3.*eosMDJWFnum(3)*t1)
247 & + eosMDJWFnum(5)*s1
248 & + p1t1*(2.*eosMDJWFnum(8) + 2.*eosMDJWFnum(11)*p1)
249
250 dden_dtheta = eosMDJWFden(1)
251 & + t1*(2.*eosMDJWFden(2)
252 & + t1*(3.*eosMDJWFden(3)
253 & + 4.*eosMDJWFden(4)*t1 ) )
254 & + s1*(eosMDJWFden(6)
255 & + t1*(3.*eosMDJWFden(7)*t1
256 & + 2.*eosMDJWFden(9)*sp5 ) )
257 & + p1*p1*(3.*eosMDJWFden(11)*t2 + eosMDJWFden(12)*p1)
258
259 alphaLoc(i,j) = rhoDen(i,j)*(dnum_dtheta
260 & - (rhoLoc(i,j)*rhoDen(i,j))*dden_dtheta)
261
262 ENDDO
263 ENDDO
264
265 ELSE
266 WRITE(*,*) 'FIND_ALPHA: equationOfState = ',equationOfState
267 STOP 'FIND_ALPHA: "equationOfState" has illegal value'
268 ENDIF
269
270 RETURN
271 END
272
273 SUBROUTINE FIND_BETA (
274 I bi, bj, iMin, iMax, jMin, jMax, k, kRef,
275 O betaloc,
276 I myThid )
277 C /==========================================================\
278 C | o SUBROUTINE FIND_BETA |
279 C | Calculates [drho(S,T,z) / dS] of a horizontal slice |
280 C |==========================================================|
281 C | |
282 C | k - is the Theta/Salt level |
283 C | kRef - determines pressure reference level |
284 C | (not used in 'LINEAR' mode) |
285 C | |
286 C | betaloc - drho / dS (kg/m^3/PSU) |
287 C | |
288 C \==========================================================/
289 IMPLICIT NONE
290
291 C === Global variables ===
292 #include "SIZE.h"
293 #include "DYNVARS.h"
294 #include "EEPARAMS.h"
295 #include "PARAMS.h"
296 #include "EOS.h"
297 #include "GRID.h"
298
299 C == Routine arguments ==
300 C k :: Level of Theta/Salt slice
301 C kRef :: Pressure reference level
302 c myThid :: thread number for this instance of the routine
303 INTEGER myThid
304 INTEGER bi,bj,iMin,iMax,jMin,jMax
305 INTEGER k
306 INTEGER kRef
307 _RL betaloc(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
308
309 C == Local variables ==
310 INTEGER i,j
311 _RL refTemp,refSalt,tP,sP
312 _RL t1, t2, t3, s1, s3o2, p1, sp5, p1t1
313 _RL drhoP0dS
314 _RL dKdS, dKdSSalt, dKdSPres
315 _RL locPres(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
316 _RL rhoP0 (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
317 _RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
318 _RL dnum_dsalt, dden_dsalt
319 _RL rhoDen (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
320 _RL rhoLoc (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
321 CEOP
322
323 #ifdef CHECK_SALINITY_FOR_NEGATIVE_VALUES
324 CALL LOOK_FOR_NEG_SALINITY( bi, bj, iMin, iMax, jMin, jMax, k,
325 & sFld, myThid )
326 #endif
327
328 IF (equationOfState.EQ.'LINEAR') THEN
329
330 DO j=jMin,jMax
331 DO i=iMin,iMax
332 betaloc(i,j) = rhonil * sBeta
333 ENDDO
334 ENDDO
335
336 ELSEIF (equationOfState.EQ.'POLY3') THEN
337
338 refTemp=eosRefT(kRef)
339 refSalt=eosRefS(kRef)
340
341 DO j=jMin,jMax
342 DO i=iMin,iMax
343 tP=theta(i,j,k,bi,bj)-refTemp
344 sP=salt(i,j,k,bi,bj)-refSalt
345 #ifdef USE_FACTORIZED_POLY
346 betaloc(i,j) =
347 & ( eosC(9,kRef)*sP*3. + eosC(5,kRef)*2. )*sP + eosC(2,kRef)
348 & + ( eosC(7,kRef)*tP
349 & +eosC(8,kRef)*sP*2. + eosC(4,kRef)
350 & )*tP
351 #else
352 betaloc(i,j) =
353 & eosC(2,kRef) +
354 & eosC(4,kRef)*tP +
355 & eosC(5,kRef) *sP*2. +
356 & eosC(7,kRef)*tP*tP +
357 & eosC(8,kRef)*tP *sP*2. +
358 & eosC(9,kRef) *sP*sP*3.
359 #endif
360 ENDDO
361 ENDDO
362
363 ELSEIF ( equationOfState(1:5).EQ.'JMD95'
364 & .OR. equationOfState.EQ.'UNESCO' ) THEN
365 C nonlinear equation of state in pressure coordinates
366
367 CALL PRESSURE_FOR_EOS(
368 I bi, bj, iMin, iMax, jMin, jMax, kRef,
369 O locPres,
370 I myThid )
371
372 CALL FIND_RHOP0(
373 I bi, bj, iMin, iMax, jMin, jMax, k,
374 I theta, salt,
375 O rhoP0,
376 I myThid )
377
378 CALL FIND_BULKMOD(
379 I bi, bj, iMin, iMax, jMin, jMax, k,
380 I locPres, theta, salt,
381 O bulkMod,
382 I myThid )
383
384 DO j=jMin,jMax
385 DO i=iMin,iMax
386
387 C abbreviations
388 t1 = theta(i,j,k,bi,bj)
389 t2 = t1*t1
390 t3 = t2*t1
391
392 s1 = salt(i,j,k,bi,bj)
393 IF ( s1 .GT. 0. _d 0 ) THEN
394 s3o2 = 1.5*SQRT(s1)
395 ELSE
396 s1 = 0. _d 0
397 s3o2 = 0. _d 0
398 ENDIF
399
400 p1 = locPres(i,j)*SItoBar
401
402 C d(rho)/d(S)
403 C of fresh water at p = 0
404 drhoP0dS = 0. _d 0
405 C of salt water at p = 0
406 drhoP0dS = drhoP0dS
407 & + eosJMDCSw(1)
408 & + eosJMDCSw(2)*t1
409 & + eosJMDCSw(3)*t2
410 & + eosJMDCSw(4)*t3
411 & + eosJMDCSw(5)*t3*t1
412 & + s3o2*(
413 & eosJMDCSw(6)
414 & + eosJMDCSw(7)*t1
415 & + eosJMDCSw(8)*t2
416 & )
417 & + 2*eosJMDCSw(9)*s1
418 C d(bulk modulus)/d(S)
419 C of fresh water at p = 0
420 dKdS = 0. _d 0
421 C of sea water at p = 0
422 dKdSSalt =
423 & eosJMDCKSw(1)
424 & + eosJMDCKSw(2)*t1
425 & + eosJMDCKSw(3)*t2
426 & + eosJMDCKSw(4)*t3
427 & + s3o2*( eosJMDCKSw(5)
428 & + eosJMDCKSw(6)*t1
429 & + eosJMDCKSw(7)*t2
430 & )
431
432 C of sea water at p
433 dKdSPres =
434 & p1*( eosJMDCKP(5)
435 & + eosJMDCKP(6)*t1
436 & + eosJMDCKP(7)*t2
437 & )
438 & + s3o2*p1*eosJMDCKP(8)
439 & + p1*p1*( eosJMDCKP(12)
440 & + eosJMDCKP(13)*t1
441 & + eosJMDCKP(14)*t2
442 & )
443
444 dKdS = dKdSSalt + dKdSPres
445
446 betaloc(i,j) =
447 & ( bulkmod(i,j)**2*drhoP0dS
448 & - bulkmod(i,j)*p1*drhoP0dS
449 & - rhoP0(i,j)*p1*dKdS )
450 & /( bulkmod(i,j) - p1 )**2
451
452
453 ENDDO
454 ENDDO
455 ELSEIF ( equationOfState.EQ.'MDJWF' ) THEN
456
457 CALL PRESSURE_FOR_EOS(
458 I bi, bj, iMin, iMax, jMin, jMax, kRef,
459 O locPres,
460 I myThid )
461
462 c$$$ CALL FIND_RHO( bi, bj, iMin, iMax, jMin, jMax, k,
463 c$$$ & kRef, theta, salt, rhoLoc, myThid )
464
465 CALL FIND_RHONUM( bi, bj, iMin, iMax, jMin, jMax, k,
466 & locPres, theta, salt, rhoLoc, myThid )
467
468 CALL FIND_RHODEN( bi, bj, iMin, iMax, jMin, jMax, k,
469 & locPres, theta, salt, rhoDen, myThid )
470
471 DO j=jMin,jMax
472 DO i=iMin,iMax
473 t1 = theta(i,j,k,bi,bj)
474 t2 = t1*t1
475 s1 = salt(i,j,k,bi,bj)
476 IF ( s1 .GT. 0. _d 0 ) THEN
477 sp5 = SQRT(s1)
478 ELSE
479 s1 = 0. _d 0
480 sp5 = 0. _d 0
481 ENDIF
482
483 p1 = locPres(i,j)*SItodBar
484 p1t1 = p1*t1
485
486 dnum_dsalt = eosMDJWFnum(4)
487 & + eosMDJWFnum(5)*t1
488 & + 2.*eosMDJWFnum(6)*s1 + eosMDJWFnum(9)*p1
489 dden_dsalt = eosMDJWFden(5)
490 & + t1*( eosMDJWFden(6) + eosMDJWFden(7)*t2 )
491 & + 1.5*sp5*(eosMDJWFden(8) + eosMDJWFden(9)*t2)
492
493 betaLoc(i,j) = rhoDen(i,j)*( dnum_dsalt
494 & - (rhoLoc(i,j)*rhoDen(i,j))*dden_dsalt )
495
496 ENDDO
497 ENDDO
498
499 ELSE
500 WRITE(*,*) 'FIND_BETA: equationOfState = ',equationOfState
501 STOP 'FIND_BETA: "equationOfState" has illegal value'
502 ENDIF
503
504 RETURN
505 END

  ViewVC Help
Powered by ViewVC 1.1.22