/[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.17 - (show annotations) (download)
Sat Aug 9 01:02:28 2008 UTC (15 years, 8 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint62v, checkpoint62u, checkpoint62t, checkpoint62c, checkpoint62s, checkpoint62r, checkpoint62q, checkpoint62p, checkpoint62a, checkpoint62g, checkpoint62f, checkpoint62e, checkpoint62d, checkpoint62k, checkpoint62j, checkpoint62i, checkpoint62h, checkpoint62o, checkpoint62n, checkpoint62m, checkpoint62l, checkpoint62w, checkpoint62z, checkpoint62y, checkpoint62x, checkpoint62, checkpoint63, checkpoint62b, checkpoint61f, checkpoint61n, checkpoint61q, checkpoint61e, checkpoint61g, checkpoint61d, checkpoint61c, checkpoint61o, checkpoint61l, checkpoint61m, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i, checkpoint61v, checkpoint61w, checkpoint61t, checkpoint61u, checkpoint61r, checkpoint61s, checkpoint61p, checkpoint61z, checkpoint61x, checkpoint61y
Changes since 1.16: +109 -99 lines
new S/R FIND_RHO_2D (replace FIND_RHO):
 change argument tFld & sFld to just 2-D local-tile slice arrays
  (+ same changes in all S/R called by FIND_RHO)
for now, keep old S/R FIND_RHO (which just call new S/R FIND_RHO_2D)

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

  ViewVC Help
Powered by ViewVC 1.1.22