/[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.15 - (show annotations) (download)
Mon Oct 27 22:32:55 2003 UTC (20 years, 7 months ago) by heimbach
Branch: MAIN
CVS Tags: checkpoint57m_post, checkpoint52l_pre, hrcube4, hrcube5, checkpoint57g_pre, checkpoint57s_post, checkpoint58b_post, checkpoint57b_post, checkpoint52d_pre, checkpoint57g_post, checkpoint56b_post, checkpoint57y_post, checkpoint52j_pre, checkpoint54d_post, checkpoint54e_post, checkpoint57r_post, checkpoint57d_post, checkpoint57i_post, checkpoint52l_post, checkpoint52k_post, checkpoint59, checkpoint58, checkpoint55, checkpoint54, checkpoint57, checkpoint56, checkpoint53, checkpoint52, checkpoint58f_post, checkpoint52f_post, checkpoint57n_post, checkpoint58d_post, checkpoint58a_post, checkpoint57z_post, checkpoint54f_post, checkpoint58y_post, checkpoint51t_post, checkpoint58t_post, checkpoint55i_post, checkpoint58m_post, checkpoint57l_post, checkpoint52i_pre, hrcube_1, hrcube_2, hrcube_3, checkpoint51s_post, checkpoint57t_post, checkpoint55c_post, checkpoint52e_pre, checkpoint57v_post, checkpoint57f_post, checkpoint52e_post, checkpoint53d_post, checkpoint57a_post, checkpoint57h_pre, checkpoint52b_pre, checkpoint54b_post, checkpoint58w_post, checkpoint57h_post, checkpoint52m_post, checkpoint57y_pre, checkpoint55g_post, checkpoint51q_post, checkpoint52b_post, checkpoint52c_post, checkpoint58o_post, checkpoint57c_post, checkpoint58p_post, checkpoint58q_post, checkpoint52f_pre, checkpoint55d_post, checkpoint58e_post, mitgcm_mapl_00, checkpoint54a_pre, checkpoint53c_post, checkpoint55d_pre, checkpoint57c_pre, checkpoint58r_post, checkpoint55j_post, checkpoint54a_post, checkpoint55h_post, checkpoint58n_post, checkpoint51r_post, checkpoint57e_post, checkpoint55b_post, checkpoint53a_post, checkpoint55f_post, checkpoint52d_post, checkpoint53g_post, checkpoint57p_post, checkpint57u_post, checkpoint57q_post, eckpoint57e_pre, checkpoint58k_post, checkpoint52a_pre, checkpoint58v_post, checkpoint52i_post, checkpoint52h_pre, checkpoint56a_post, checkpoint58l_post, checkpoint53f_post, checkpoint57h_done, checkpoint52j_post, checkpoint57j_post, checkpoint57f_pre, checkpoint58g_post, branch-netcdf, checkpoint58x_post, checkpoint52n_post, checkpoint53b_pre, checkpoint58h_post, checkpoint56c_post, checkpoint58j_post, checkpoint57a_pre, checkpoint55a_post, checkpoint57o_post, checkpoint51o_post, checkpoint57k_post, checkpoint53b_post, checkpoint52a_post, checkpoint57w_post, checkpoint58i_post, ecco_c52_e35, checkpoint57x_post, checkpoint58c_post, checkpoint58u_post, checkpoint53d_pre, checkpoint58s_post, checkpoint55e_post, checkpoint54c_post, checkpoint51p_post, checkpoint51u_post
Branch point for: branch-nonh, netcdf-sm0
Changes since 1.14: +13 -7 lines
o cleaning ALLOW_GRADIENT_CHECK -> ALLOW_GRDCHK
o cleaning some ALLOW_TANGENTLINEAR_RUN -> ALLOW_AUTODIFF
o bug fix in find_alpha.F for MDJWF:
  - modif. to alpha = 1/D*( dN/dT - rho*dD/Dt) to account for
    change rho -> rho-rhoConst
  - replace call find_rho to find_rhonum

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

  ViewVC Help
Powered by ViewVC 1.1.22