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

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

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


Revision 1.26 - (show annotations) (download)
Mon Sep 29 19:24:31 2003 UTC (20 years, 7 months ago) by edhill
Branch: MAIN
CVS Tags: checkpoint51k_post, checkpoint52l_pre, hrcube4, hrcube5, checkpoint57b_post, checkpoint52d_pre, checkpoint56b_post, checkpoint52j_pre, checkpoint51o_pre, checkpoint54d_post, checkpoint54e_post, checkpoint51l_post, checkpoint52l_post, checkpoint52k_post, checkpoint55, checkpoint54, checkpoint57, checkpoint56, checkpoint53, checkpoint52, checkpoint52f_post, checkpoint54f_post, checkpoint51f_post, checkpoint51t_post, checkpoint51n_post, checkpoint55i_post, checkpoint52i_pre, hrcube_1, hrcube_2, hrcube_3, checkpoint51s_post, checkpoint55c_post, checkpoint51j_post, checkpoint52e_pre, checkpoint52e_post, checkpoint51n_pre, checkpoint53d_post, checkpoint57a_post, checkpoint52b_pre, checkpoint54b_post, checkpoint51l_pre, checkpoint52m_post, checkpoint55g_post, checkpoint51q_post, checkpoint52b_post, checkpoint52c_post, checkpoint51h_pre, checkpoint52f_pre, checkpoint55d_post, checkpoint54a_pre, checkpoint53c_post, checkpoint55d_pre, checkpoint55j_post, checkpoint54a_post, checkpoint55h_post, checkpoint51r_post, checkpoint51i_post, checkpoint55b_post, checkpoint53a_post, checkpoint55f_post, checkpoint52d_post, checkpoint53g_post, checkpoint52a_pre, checkpoint52i_post, checkpoint51i_pre, checkpoint52h_pre, checkpoint56a_post, checkpoint53f_post, checkpoint52j_post, branch-netcdf, checkpoint52n_post, checkpoint53b_pre, checkpoint56c_post, checkpoint57a_pre, checkpoint55a_post, checkpoint51o_post, checkpoint53b_post, checkpoint52a_post, checkpoint51g_post, ecco_c52_e35, checkpoint51m_post, checkpoint53d_pre, checkpoint55e_post, checkpoint54c_post, checkpoint51p_post, checkpoint51u_post
Branch point for: branch-nonh, tg2-branch, netcdf-sm0, checkpoint51n_branch
Changes since 1.25: +3 -3 lines
 o convert all comments with single quotes (such as "can't", "don't", etc.)
     to unabbreviated form (eg. "do not") since these unmatched quotes
     tend to break the Fortran parser used to generate the HTML-ified
     code browser (see: http://mitgcm.org/dev_docs/code_reference/)

1 C $Header: /u/u3/gcmpack/MITgcm/model/src/find_rho.F,v 1.25 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_RHO
9 C !INTERFACE:
10 SUBROUTINE FIND_RHO(
11 I bi, bj, iMin, iMax, jMin, jMax, k, kRef,
12 I tFld, sFld,
13 O rholoc,
14 I myThid )
15
16 C !DESCRIPTION: \bv
17 C *==========================================================*
18 C | o SUBROUTINE FIND_RHO
19 C | Calculates [rho(S,T,z)-rhoConst] of a slice
20 C *==========================================================*
21 C |
22 C | k - is the Theta/Salt level
23 C | kRef - determines pressure reference level
24 C | (not used in 'LINEAR' mode)
25 C |
26 C *==========================================================*
27 C \ev
28
29 C !USES:
30 IMPLICIT NONE
31 C == Global variables ==
32 #include "SIZE.h"
33 #include "EEPARAMS.h"
34 #include "PARAMS.h"
35 #include "EOS.h"
36 #include "GRID.h"
37
38 C !INPUT/OUTPUT PARAMETERS:
39 C == Routine arguments ==
40 C k :: Level of Theta/Salt slice
41 C kRef :: Pressure reference level
42 INTEGER bi,bj,iMin,iMax,jMin,jMax
43 INTEGER k
44 INTEGER kRef
45 _RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
46 _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
47 _RL rholoc(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
48 INTEGER myThid
49
50 C !LOCAL VARIABLES:
51 C == Local variables ==
52 INTEGER i,j
53 _RL refTemp,refSalt,sigRef,tP,sP,deltaSig,dRho
54 _RL locPres(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
55 _RL rhoP0 (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
56 _RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
57 _RL rhoNum (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
58 _RL rhoDen (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
59 CHARACTER*(MAX_LEN_MBUF) msgbuf
60 CEOP
61
62 #ifdef ALLOW_AUTODIFF_TAMC
63 DO j=1-OLy,sNy+OLy
64 DO i=1-OLx,sNx+OLx
65 rholoc(i,j) = 0. _d 0
66 rhoP0(i,j) = 0. _d 0
67 bulkMod(i,j) = 0. _d 0
68 ENDDO
69 ENDDO
70 #endif
71
72 #ifdef CHECK_SALINITY_FOR_NEGATIVE_VALUES
73 CALL LOOK_FOR_NEG_SALINITY( bi, bj, iMin, iMax, jMin, jMax, k,
74 & sFld, myThid )
75 #endif
76
77 IF (equationOfState.EQ.'LINEAR') THEN
78
79 C ***NOTE***
80 C In the linear EOS, to make the static stability calculation meaningful
81 C we alway calculate the perturbation with respect to the surface level.
82 C **********
83 refTemp=tRef(kRef)
84 refSalt=sRef(kRef)
85
86 dRho = rhoNil-rhoConst
87
88 DO j=jMin,jMax
89 DO i=iMin,iMax
90 rholoc(i,j)=rhoNil*(
91 & sBeta*(sFld(i,j,k,bi,bj)-refSalt)
92 & -tAlpha*(tFld(i,j,k,bi,bj)-refTemp) )
93 & + dRho
94 ENDDO
95 ENDDO
96
97 ELSEIF (equationOfState.EQ.'POLY3') THEN
98
99 refTemp=eosRefT(kRef)
100 refSalt=eosRefS(kRef)
101 sigRef=eosSig0(kRef) + (1000.-rhoConst)
102
103 DO j=jMin,jMax
104 DO i=iMin,iMax
105 tP=tFld(i,j,k,bi,bj)-refTemp
106 sP=sFld(i,j,k,bi,bj)-refSalt
107 #ifdef USE_FACTORIZED_POLY
108 deltaSig=
109 & (( eosC(9,kRef)*sP + eosC(5,kRef) )*sP + eosC(2,kRef) )*sP
110 & + ( ( eosC(6,kRef)
111 & *tP
112 & +eosC(7,kRef)*sP + eosC(3,kRef)
113 & )*tP
114 & +(eosC(8,kRef)*sP + eosC(4,kRef) )*sP + eosC(1,kRef)
115 & )*tP
116 #else
117 deltaSig=
118 & eosC(1,kRef)*tP
119 & +eosC(2,kRef) *sP
120 & +eosC(3,kRef)*tP*tP
121 & +eosC(4,kRef)*tP *sP
122 & +eosC(5,kRef) *sP*sP
123 & +eosC(6,kRef)*tP*tP*tP
124 & +eosC(7,kRef)*tP*tP *sP
125 & +eosC(8,kRef)*tP *sP*sP
126 & +eosC(9,kRef) *sP*sP*sP
127 #endif
128 rholoc(i,j)=sigRef+deltaSig
129 ENDDO
130 ENDDO
131
132 ELSEIF ( equationOfState(1:5).EQ.'JMD95'
133 & .OR. equationOfState.EQ.'UNESCO' ) THEN
134 C nonlinear equation of state in pressure coordinates
135
136 CALL PRESSURE_FOR_EOS(
137 I bi, bj, iMin, iMax, jMin, jMax, kRef,
138 O locPres,
139 I myThid )
140
141 CALL FIND_RHOP0(
142 I bi, bj, iMin, iMax, jMin, jMax, k,
143 I tFld, sFld,
144 O rhoP0,
145 I myThid )
146
147 CALL FIND_BULKMOD(
148 I bi, bj, iMin, iMax, jMin, jMax, k,
149 I locPres, tFld, sFld,
150 O bulkMod,
151 I myThid )
152
153 #ifdef ALLOW_AUTODIFF_TAMC
154 cph can not DO storing here since find_rho is called multiple times;
155 cph additional recomp. should be acceptable
156 cphCADJ STORE rhoP0(:,:) = comlev1_bibj_k , key=kkey , byte=isbyte
157 cphCADJ STORE bulkMod(:,:) = comlev1_bibj_k , key=kkey , byte=isbyte
158 #endif /* ALLOW_AUTODIFF_TAMC */
159 DO j=jMin,jMax
160 DO i=iMin,iMax
161
162 C density of sea water at pressure p
163 rholoc(i,j) = rhoP0(i,j)
164 & /(1. _d 0 -
165 & locPres(i,j)*SItoBar/bulkMod(i,j))
166 & - rhoConst
167
168 ENDDO
169 ENDDO
170
171 ELSEIF ( equationOfState.EQ.'MDJWF' ) THEN
172
173 CALL PRESSURE_FOR_EOS(
174 I bi, bj, iMin, iMax, jMin, jMax, kRef,
175 O locPres,
176 I myThid )
177
178 CALL FIND_RHONUM( bi, bj, iMin, iMax, jMin, jMax, k,
179 & locPres, tFld, sFld, rhoNum, myThid )
180
181 CALL FIND_RHODEN( bi, bj, iMin, iMax, jMin, jMax, k,
182 & locPres, tFld, sFld, rhoDen, myThid )
183 #ifdef ALLOW_AUTODIFF_TAMC
184 cph can not DO storing here since find_rho is called multiple times;
185 cph additional recomp. should be acceptable
186 cphCADJ STORE rhoNum(:,:) = comlev1_bibj_k , key=kkey , byte=isbyte
187 cphCADJ STORE rhoDen(:,:) = comlev1_bibj_k , key=kkey , byte=isbyte
188 #endif /* ALLOW_AUTODIFF_TAMC */
189 DO j=jMin,jMax
190 DO i=iMin,iMax
191 rholoc(i,j) = rhoNum(i,j)*rhoDen(i,j) - rhoConst
192 ENDDO
193 ENDDO
194
195 ELSEIF( equationOfState .EQ. 'IDEALG' ) THEN
196 C
197 ELSE
198 WRITE(msgbuf,'(3a)')
199 & ' FIND_RHO: equationOfState = "',equationOfState,'"'
200 CALL print_error( msgbuf, mythid )
201 STOP 'ABNORMAL END: S/R FIND_RHO'
202 ENDIF
203
204 RETURN
205 END
206
207 CBOP
208 C !ROUTINE: FIND_RHOP0
209 C !INTERFACE:
210 SUBROUTINE FIND_RHOP0(
211 I bi, bj, iMin, iMax, jMin, jMax, k,
212 I tFld, sFld,
213 O rhoP0,
214 I myThid )
215
216 C !DESCRIPTION: \bv
217 C *==========================================================*
218 C | o SUBROUTINE FIND_RHOP0
219 C | Calculates rho(S,T,0) of a slice
220 C *==========================================================*
221 C |
222 C | k - is the surface level
223 C |
224 C *==========================================================*
225 C \ev
226
227 C !USES:
228 IMPLICIT NONE
229 C == Global variables ==
230 #include "SIZE.h"
231 #include "EEPARAMS.h"
232 #include "PARAMS.h"
233 #include "EOS.h"
234
235 C !INPUT/OUTPUT PARAMETERS:
236 C == Routine arguments ==
237 C k :: Level of Theta/Salt slice
238 INTEGER bi,bj,iMin,iMax,jMin,jMax
239 INTEGER k
240 _RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
241 _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
242 _RL rhoP0(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
243 _RL rfresh, rsalt
244 INTEGER myThid
245
246 C !LOCAL VARIABLES:
247 C == Local variables ==
248 INTEGER i,j
249 _RL t, t2, t3, t4, s, s3o2
250 CEOP
251
252 DO j=jMin,jMax
253 DO i=iMin,iMax
254 C abbreviations
255 t = tFld(i,j,k,bi,bj)
256 t2 = t*t
257 t3 = t2*t
258 t4 = t3*t
259
260 s = sFld(i,j,k,bi,bj)
261 IF ( s .GT. 0. _d 0 ) THEN
262 s3o2 = s*SQRT(s)
263 ELSE
264 s = 0. _d 0
265 s3o2 = 0. _d 0
266 ENDIF
267
268 C density of freshwater at the surface
269 rfresh =
270 & eosJMDCFw(1)
271 & + eosJMDCFw(2)*t
272 & + eosJMDCFw(3)*t2
273 & + eosJMDCFw(4)*t3
274 & + eosJMDCFw(5)*t4
275 & + eosJMDCFw(6)*t4*t
276 C density of sea water at the surface
277 rsalt =
278 & s*(
279 & eosJMDCSw(1)
280 & + eosJMDCSw(2)*t
281 & + eosJMDCSw(3)*t2
282 & + eosJMDCSw(4)*t3
283 & + eosJMDCSw(5)*t4
284 & )
285 & + s3o2*(
286 & eosJMDCSw(6)
287 & + eosJMDCSw(7)*t
288 & + eosJMDCSw(8)*t2
289 & )
290 & + eosJMDCSw(9)*s*s
291
292 rhoP0(i,j) = rfresh + rsalt
293
294 ENDDO
295 ENDDO
296
297 RETURN
298 END
299
300 C !ROUTINE: FIND_BULKMOD
301 C !INTERFACE:
302 SUBROUTINE FIND_BULKMOD(
303 I bi, bj, iMin, iMax, jMin, jMax, k,
304 I locPres, tFld, sFld,
305 O bulkMod,
306 I myThid )
307
308 C !DESCRIPTION: \bv
309 C *==========================================================*
310 C | o SUBROUTINE FIND_BULKMOD
311 C | Calculates the secant bulk modulus K(S,T,p) of a slice
312 C *==========================================================*
313 C |
314 C | k - is the level of Theta/Salt slice
315 C |
316 C *==========================================================*
317 C \ev
318
319 C !USES:
320 IMPLICIT NONE
321 C == Global variables ==
322 #include "SIZE.h"
323 #include "EEPARAMS.h"
324 #include "PARAMS.h"
325 #include "EOS.h"
326
327 C !INPUT/OUTPUT PARAMETERS:
328 C == Routine arguments ==
329 C k :: Level of Theta/Salt slice
330 INTEGER bi,bj,iMin,iMax,jMin,jMax
331 INTEGER k
332 INTEGER kRef
333 _RL locPres(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
334 _RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
335 _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
336 _RL bulkMod(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
337 _RL bMfresh, bMsalt, bMpres
338 INTEGER myThid
339
340 C !LOCAL VARIABLES:
341 C == Local variables ==
342 INTEGER i,j
343 _RL t, t2, t3, t4, s, s3o2, p, p2
344 CEOP
345
346 DO j=jMin,jMax
347 DO i=iMin,iMax
348 C abbreviations
349 t = tFld(i,j,k,bi,bj)
350 t2 = t*t
351 t3 = t2*t
352 t4 = t3*t
353
354 s = sFld(i,j,k,bi,bj)
355 IF ( s .GT. 0. _d 0 ) THEN
356 s3o2 = s*SQRT(s)
357 ELSE
358 s = 0. _d 0
359 s3o2 = 0. _d 0
360 ENDIF
361 C
362 p = locPres(i,j)*SItoBar
363 p2 = p*p
364 C secant bulk modulus of fresh water at the surface
365 bMfresh =
366 & eosJMDCKFw(1)
367 & + eosJMDCKFw(2)*t
368 & + eosJMDCKFw(3)*t2
369 & + eosJMDCKFw(4)*t3
370 & + eosJMDCKFw(5)*t4
371 C secant bulk modulus of sea water at the surface
372 bMsalt =
373 & s*( eosJMDCKSw(1)
374 & + eosJMDCKSw(2)*t
375 & + eosJMDCKSw(3)*t2
376 & + eosJMDCKSw(4)*t3
377 & )
378 & + s3o2*( eosJMDCKSw(5)
379 & + eosJMDCKSw(6)*t
380 & + eosJMDCKSw(7)*t2
381 & )
382 C secant bulk modulus of sea water at pressure p
383 bMpres =
384 & p*( eosJMDCKP(1)
385 & + eosJMDCKP(2)*t
386 & + eosJMDCKP(3)*t2
387 & + eosJMDCKP(4)*t3
388 & )
389 & + p*s*( eosJMDCKP(5)
390 & + eosJMDCKP(6)*t
391 & + eosJMDCKP(7)*t2
392 & )
393 & + p*s3o2*eosJMDCKP(8)
394 & + p2*( eosJMDCKP(9)
395 & + eosJMDCKP(10)*t
396 & + eosJMDCKP(11)*t2
397 & )
398 & + p2*s*( eosJMDCKP(12)
399 & + eosJMDCKP(13)*t
400 & + eosJMDCKP(14)*t2
401 & )
402
403 bulkMod(i,j) = bMfresh + bMsalt + bMpres
404
405 ENDDO
406 ENDDO
407
408 RETURN
409 END
410
411 CBOP
412 C !ROUTINE: FIND_RHONUM
413 C !INTERFACE:
414 SUBROUTINE FIND_RHONUM(
415 I bi, bj, iMin, iMax, jMin, jMax, k,
416 I locPres, tFld, sFld,
417 O rhoNum,
418 I myThid )
419
420 C !DESCRIPTION: \bv
421 C *==========================================================*
422 C | o SUBROUTINE FIND_RHONUM
423 C | Calculates the numerator of the McDougall et al.
424 C | equation of state
425 C | - the code is more or less a copy of MOM4
426 C *==========================================================*
427 C |
428 C | k - is the level of Theta/Salt slice
429 C |
430 C *==========================================================*
431 C \ev
432
433 C !USES:
434 IMPLICIT NONE
435 C == Global variables ==
436 #include "SIZE.h"
437 #include "EEPARAMS.h"
438 #include "PARAMS.h"
439 #include "EOS.h"
440
441 C !INPUT/OUTPUT PARAMETERS:
442 C == Routine arguments ==
443 C k :: Level of Theta/Salt slice
444 INTEGER bi,bj,iMin,iMax,jMin,jMax
445 INTEGER k
446 _RL locPres(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
447 _RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
448 _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
449 _RL rhoNum(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
450 INTEGER myThid
451
452 C !LOCAL VARIABLES:
453 C == Local variables ==
454 INTEGER i,j
455 _RL t1, t2, s1, p1
456 CEOP
457 DO j=jMin,jMax
458 DO i=iMin,iMax
459 C abbreviations
460 t1 = tFld(i,j,k,bi,bj)
461 t2 = t1*t1
462 s1 = sFld(i,j,k,bi,bj)
463
464 p1 = locPres(i,j)*SItodBar
465
466 rhoNum(i,j) = eosMDJWFnum(0)
467 & + t1*(eosMDJWFnum(1)
468 & + t1*(eosMDJWFnum(2) + eosMDJWFnum(3)*t1) )
469 & + s1*(eosMDJWFnum(4)
470 & + eosMDJWFnum(5)*t1 + eosMDJWFnum(6)*s1)
471 & + p1*(eosMDJWFnum(7) + eosMDJWFnum(8)*t2
472 & + eosMDJWFnum(9)*s1
473 & + p1*(eosMDJWFnum(10) + eosMDJWFnum(11)*t2) )
474
475 ENDDO
476 ENDDO
477
478 RETURN
479 end
480
481 CBOP
482 C !ROUTINE: FIND_RHODEN
483 C !INTERFACE:
484 SUBROUTINE FIND_RHODEN(
485 I bi, bj, iMin, iMax, jMin, jMax, k,
486 I locPres, tFld, sFld,
487 O rhoDen,
488 I myThid )
489
490 C !DESCRIPTION: \bv
491 C *==========================================================*
492 C | o SUBROUTINE FIND_RHODEN
493 C | Calculates the denominator of the McDougall et al.
494 C | equation of state
495 C | - the code is more or less a copy of MOM4
496 C *==========================================================*
497 C |
498 C | k - is the level of Theta/Salt slice
499 C |
500 C *==========================================================*
501 C \ev
502
503 C !USES:
504 IMPLICIT NONE
505 C == Global variables ==
506 #include "SIZE.h"
507 #include "EEPARAMS.h"
508 #include "PARAMS.h"
509 #include "EOS.h"
510
511 C !INPUT/OUTPUT PARAMETERS:
512 C == Routine arguments ==
513 C k :: Level of Theta/Salt slice
514 INTEGER bi,bj,iMin,iMax,jMin,jMax
515 INTEGER k
516 _RL locPres(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
517 _RL tFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
518 _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
519 _RL rhoDen(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
520 INTEGER myThid
521
522 C !LOCAL VARIABLES:
523 C == Local variables ==
524 INTEGER i,j
525 _RL t1, t2, s1, sp5, p1, p1t1
526 _RL den, epsln
527 parameter ( epsln = 0. _d 0 )
528 CEOP
529 DO j=jMin,jMax
530 DO i=iMin,iMax
531 C abbreviations
532 t1 = tFld(i,j,k,bi,bj)
533 t2 = t1*t1
534 s1 = sFld(i,j,k,bi,bj)
535 IF ( s1 .GT. 0. _d 0 ) THEN
536 sp5 = SQRT(s1)
537 ELSE
538 s1 = 0. _d 0
539 sp5 = 0. _d 0
540 ENDIF
541
542 p1 = locPres(i,j)*SItodBar
543 p1t1 = p1*t1
544
545 den = eosMDJWFden(0)
546 & + t1*(eosMDJWFden(1)
547 & + t1*(eosMDJWFden(2)
548 & + t1*(eosMDJWFden(3) + t1*eosMDJWFden(4) ) ) )
549 & + s1*(eosMDJWFden(5)
550 & + t1*(eosMDJWFden(6)
551 & + eosMDJWFden(7)*t2)
552 & + sp5*(eosMDJWFden(8) + eosMDJWFden(9)*t2) )
553 & + p1*(eosMDJWFden(10)
554 & + p1t1*(eosMDJWFden(11)*t2 + eosMDJWFden(12)*p1) )
555
556 rhoDen(i,j) = 1.0/(epsln+den)
557
558 ENDDO
559 ENDDO
560
561 RETURN
562 end
563
564 SUBROUTINE find_rho_scalar(
565 I tLoc, sLoc, pLoc,
566 O rhoLoc,
567 I myThid )
568
569 C !DESCRIPTION: \bv
570 C *==========================================================*
571 C | o SUBROUTINE FIND_RHO_SCALAR
572 C | Calculates [rho(S,T,p)-rhoConst]
573 C *==========================================================*
574 C \ev
575
576 C !USES:
577 IMPLICIT NONE
578 C == Global variables ==
579 #include "SIZE.h"
580 #include "EEPARAMS.h"
581 #include "PARAMS.h"
582 #include "EOS.h"
583
584 C !INPUT/OUTPUT PARAMETERS:
585 C == Routine arguments ==
586 _RL sLoc, tLoc, pLoc
587 _RL rhoLoc
588 INTEGER myThid
589
590 C !LOCAL VARIABLES:
591 C == Local variables ==
592
593 _RL t1, t2, t3, t4, s1, s3o2, p1, p2, sp5, p1t1
594 _RL rfresh, rsalt, rhoP0
595 _RL bMfresh, bMsalt, bMpres, BulkMod
596 _RL rhoNum, rhoDen, den, epsln
597 parameter ( epsln = 0. _d 0 )
598
599 character*(max_len_mbuf) msgbuf
600 CEOP
601
602 rhoLoc = 0. _d 0
603 rhoP0 = 0. _d 0
604 bulkMod = 0. _d 0
605 rfresh = 0. _d 0
606 rsalt = 0. _d 0
607 bMfresh = 0. _d 0
608 bMsalt = 0. _d 0
609 bMpres = 0. _d 0
610 rhoNum = 0. _d 0
611 rhoDen = 0. _d 0
612 den = 0. _d 0
613
614 t1 = tLoc
615 t2 = t1*t1
616 t3 = t2*t1
617 t4 = t3*t1
618
619 s1 = sLoc
620 IF ( s1 .LT. 0. _d 0 ) THEN
621 C issue a warning
622 WRITE(*,'(a,i3,a,i3,a,i3,a,e13.5)')
623 & ' FIND_RHO_SCALAR: WARNING, salinity = ', s1
624 s1 = 0. _d 0
625 ENDIF
626
627 IF (equationOfState.EQ.'LINEAR') THEN
628
629 rhoLoc = 0. _d 0
630
631 ELSEIF (equationOfState.EQ.'POLY3') THEN
632
633 C this is not correct, there is a field eosSig0 which should be use here
634 C but I DO not intent to include the reference level in this routine
635 WRITE(*,'(a)')
636 & ' FIND_RHO_SCALAR: for POLY3, the density is not'
637 WRITE(*,'(a)')
638 & ' computed correctly in this routine'
639 rhoLoc = 0. _d 0
640
641 ELSEIF ( equationOfState(1:5).EQ.'JMD95'
642 & .OR. equationOfState.EQ.'UNESCO' ) THEN
643 C nonlinear equation of state in pressure coordinates
644
645 s3o2 = s1*SQRT(s1)
646
647 p1 = pLoc*SItoBar
648 p2 = p1*p1
649
650 C density of freshwater at the surface
651 rfresh =
652 & eosJMDCFw(1)
653 & + eosJMDCFw(2)*t1
654 & + eosJMDCFw(3)*t2
655 & + eosJMDCFw(4)*t3
656 & + eosJMDCFw(5)*t4
657 & + eosJMDCFw(6)*t4*t1
658 C density of sea water at the surface
659 rsalt =
660 & s1*(
661 & eosJMDCSw(1)
662 & + eosJMDCSw(2)*t1
663 & + eosJMDCSw(3)*t2
664 & + eosJMDCSw(4)*t3
665 & + eosJMDCSw(5)*t4
666 & )
667 & + s3o2*(
668 & eosJMDCSw(6)
669 & + eosJMDCSw(7)*t1
670 & + eosJMDCSw(8)*t2
671 & )
672 & + eosJMDCSw(9)*s1*s1
673
674 rhoP0 = rfresh + rsalt
675
676 C secant bulk modulus of fresh water at the surface
677 bMfresh =
678 & eosJMDCKFw(1)
679 & + eosJMDCKFw(2)*t1
680 & + eosJMDCKFw(3)*t2
681 & + eosJMDCKFw(4)*t3
682 & + eosJMDCKFw(5)*t4
683 C secant bulk modulus of sea water at the surface
684 bMsalt =
685 & s1*( eosJMDCKSw(1)
686 & + eosJMDCKSw(2)*t1
687 & + eosJMDCKSw(3)*t2
688 & + eosJMDCKSw(4)*t3
689 & )
690 & + s3o2*( eosJMDCKSw(5)
691 & + eosJMDCKSw(6)*t1
692 & + eosJMDCKSw(7)*t2
693 & )
694 C secant bulk modulus of sea water at pressure p
695 bMpres =
696 & p1*( eosJMDCKP(1)
697 & + eosJMDCKP(2)*t1
698 & + eosJMDCKP(3)*t2
699 & + eosJMDCKP(4)*t3
700 & )
701 & + p1*s1*( eosJMDCKP(5)
702 & + eosJMDCKP(6)*t1
703 & + eosJMDCKP(7)*t2
704 & )
705 & + p1*s3o2*eosJMDCKP(8)
706 & + p2*( eosJMDCKP(9)
707 & + eosJMDCKP(10)*t1
708 & + eosJMDCKP(11)*t2
709 & )
710 & + p2*s1*( eosJMDCKP(12)
711 & + eosJMDCKP(13)*t1
712 & + eosJMDCKP(14)*t2
713 & )
714
715 bulkMod = bMfresh + bMsalt + bMpres
716
717 C density of sea water at pressure p
718 rhoLoc = rhoP0/(1. _d 0 - p1/bulkMod) - rhoConst
719
720 ELSEIF ( equationOfState.EQ.'MDJWF' ) THEN
721
722 sp5 = SQRT(s1)
723
724 p1 = pLoc*SItodBar
725 p1t1 = p1*t1
726
727 rhoNum = eosMDJWFnum(0)
728 & + t1*(eosMDJWFnum(1)
729 & + t1*(eosMDJWFnum(2) + eosMDJWFnum(3)*t1) )
730 & + s1*(eosMDJWFnum(4)
731 & + eosMDJWFnum(5)*t1 + eosMDJWFnum(6)*s1)
732 & + p1*(eosMDJWFnum(7) + eosMDJWFnum(8)*t2
733 & + eosMDJWFnum(9)*s1
734 & + p1*(eosMDJWFnum(10) + eosMDJWFnum(11)*t2) )
735
736
737 den = eosMDJWFden(0)
738 & + t1*(eosMDJWFden(1)
739 & + t1*(eosMDJWFden(2)
740 & + t1*(eosMDJWFden(3) + t1*eosMDJWFden(4) ) ) )
741 & + s1*(eosMDJWFden(5)
742 & + t1*(eosMDJWFden(6)
743 & + eosMDJWFden(7)*t2)
744 & + sp5*(eosMDJWFden(8) + eosMDJWFden(9)*t2) )
745 & + p1*(eosMDJWFden(10)
746 & + p1t1*(eosMDJWFden(11)*t2 + eosMDJWFden(12)*p1) )
747
748 rhoDen = 1.0/(epsln+den)
749
750 rhoLoc = rhoNum*rhoDen - rhoConst
751
752 ELSEIF( equationOfState .EQ. 'IDEALG' ) THEN
753 C
754 ELSE
755 WRITE(msgbuf,'(3A)')
756 & ' FIND_RHO_SCALAR : equationOfState = "',
757 & equationOfState,'"'
758 CALL PRINT_ERROR( msgbuf, mythid )
759 STOP 'ABNORMAL END: S/R FIND_RHO_SCALAR'
760 ENDIF
761
762 RETURN
763 END
764
765 CBOP
766 C !ROUTINE: LOOK_FOR_NEG_SALINITY
767 C !INTERFACE:
768 SUBROUTINE LOOK_FOR_NEG_SALINITY(
769 I bi, bj, iMin, iMax, jMin, jMax, k,
770 I sFld,
771 I myThid )
772
773 C !DESCRIPTION: \bv
774 C *==========================================================*
775 C | o SUBROUTINE LOOK_FOR_NEG_SALINITY
776 C | looks for and fixes negative salinity values
777 C | this is necessary IF the equation of state uses
778 C | the square root of salinity
779 C *==========================================================*
780 C |
781 C | k - is the Salt level
782 C |
783 C *==========================================================*
784 C \ev
785
786 C !USES:
787 IMPLICIT NONE
788 C == Global variables ==
789 #include "SIZE.h"
790 #include "EEPARAMS.h"
791 #include "PARAMS.h"
792 #include "GRID.h"
793
794 C !INPUT/OUTPUT PARAMETERS:
795 C == Routine arguments ==
796 C k :: Level of Theta/Salt slice
797 INTEGER bi,bj,iMin,iMax,jMin,jMax
798 INTEGER k
799 _RL sFld(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr,nSx,nSy)
800 INTEGER myThid
801
802 C !LOCAL VARIABLES:
803 C == Local variables ==
804 INTEGER i,j, localWarning
805 character*(max_len_mbuf) msgbuf
806 CEOP
807
808 localWarning = 0
809 DO j=jMin,jMax
810 DO i=iMin,iMax
811 C abbreviations
812 IF ( sFld(i,j,k,bi,bj) .LT. 0. _d 0 ) THEN
813 localWarning = localWarning + 1
814 sFld(i,j,k,bi,bj) = 0. _d 0
815 ENDIF
816 ENDDO
817 ENDDO
818 C issue a warning
819 IF ( localWarning .GT. 0 ) THEN
820 WRITE(standardMessageUnit,'(A,A)')
821 & 'S/R LOOK_FOR_NEG_SALINITY: found negative salinity',
822 & 'values and reset them to zero.'
823 WRITE(standardMessageUnit,'(A,I3)')
824 & 'S/R LOOK_FOR_NEG_SALINITY: current level is k = ', k
825 ENDIF
826
827 RETURN
828 END

  ViewVC Help
Powered by ViewVC 1.1.22